This module holds 4 types of functions:
(1) algebraic functions: "percival_seaton_coefficient",
(2) geometric functions: "triangle_inequality_holds", "is_sum_even", "zero_projections_3j_condition"
(3) bessel functions: "riccati_bessel_j", "riccati_bessel_y" and "modified_bessel_k_ratio" that call functions from special_functions.f90 library
(4) interpolation procedures: "spline" and "ispline" functions
checks if the sum of 3 integers is an even integer
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | x |
variables to check if the sum is even |
||
integer(kind=int32), | intent(in) | :: | y |
variables to check if the sum is even |
||
integer(kind=int32), | intent(in) | :: | z |
variables to check if the sum is even |
(out) result: true/false
returns interpolated value at guven u_ point number of points and ascending order of x is not checked since ispline is called after "spline" where these checks are done
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | u_ |
point at which the tabulated value is interpolated |
||
integer(kind=int32), | intent(in) | :: | N_ |
number of grid points |
||
real(kind=dp), | intent(in) | :: | x_(N_) |
grid points |
||
real(kind=dp), | intent(in) | :: | y_(N_) |
tabulated values |
||
real(kind=dp), | intent(in) | :: | b_(N_) |
arrays with coefficients of the spline function |
||
real(kind=dp), | intent(in) | :: | c_(N_) |
arrays with coefficients of the spline function |
||
real(kind=dp), | intent(in) | :: | d_(N_) |
arrays with coefficients of the spline function |
interpolated value at u_
calculates Percival-Seaton coefficients (body-fixed variant)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | j_ |
pre-collisional rotational angular momentum |
||
integer(kind=int32), | intent(in) | :: | j_prime_ |
post-collisional rotational angular momentum |
||
integer(kind=int32), | intent(in) | :: | lambda_ |
Legendre expansion coefficient \( \lambda\) |
||
integer(kind=int32), | intent(in) | :: | omega_ |
\(\bar{\Omega}\) |
(out) result: percival seaton coefficient in the body-fixed frame
Calculates 1/Gamma(1-X); modernized version of Molscat's rgamma function; see: https://github.com/molscat/molscat/blob/36fa8f93a92f851e9d84245dd6a972e2910541c5/source_code/rbesjy.f --------------------------------------------------------------------!
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | x | |||
real(kind=dp), | intent(out) | :: | odd | |||
real(kind=dp), | intent(out) | :: | even |
check if the triangle inequality for 3 variables hols
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | x |
variables to check the triangle inequality |
||
integer(kind=int32), | intent(in) | :: | y |
variables to check the triangle inequality |
||
integer(kind=int32), | intent(in) | :: | z |
variables to check the triangle inequality |
(out) result: true/false
checks the condition for nonvanishing 3-j symbol with zero projections: triangle inequality on x,y,z and if the sum x+y+z is an even integer
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | x |
variables to check for 3-j symbol conditions |
||
integer(kind=int32), | intent(in) | :: | y |
variables to check for 3-j symbol conditions |
||
integer(kind=int32), | intent(in) | :: | z |
variables to check for 3-j symbol conditions |
(out) result: true/false if conditions are met
calculates the ratio of the modified Bessel function of the second kind K_{l_ + 1/2}(x) and its first derivative (Eq. 8 in the "Solution of the coupled equations" section) Uses Temme's algorithm [N. M. Temme, J. Comput. Phys. 19 (1975) 324], implemented in "modified_bessel_temme_algorithm" subroutine; Unfortunately, the "ikv" function from special_functions library failed at large x_ values.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | l_ |
l - order of the function (without the 1/2 factor!) |
||
real(kind=dp), | intent(in) | :: | x_ |
x - argument of the function |
||
real(kind=dp), | intent(inout) | :: | ratio_ |
ratio of the modified Bessel function of the second kind to its derivative |
Implementation of the Temme's algorithm [N. M. Temme, J. Comput. Phys. 19 (1975) 324] to calculating modified Bessel functions of the second kind. This is a direct modernization of the "mbessk" subroutine in MOLSCAT: https://github.com/molscat/molscat/blob/master/source_code/rbessk.f
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | v | |||
real(kind=dp), | intent(in) | :: | x | |||
real(kind=dp), | intent(out) | :: | ck | |||
real(kind=dp), | intent(out) | :: | dk | |||
real(kind=dp), | intent(out) | :: | ek |
calculates the Riccati-Bessel function of the first kind and its first derivative. Calls the rctj function from special_functions.f90
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | l_ |
l - order of the Riccati-Bessel function of the first kind |
||
real(kind=dp), | intent(in) | :: | x_ |
x - argument of the Riccati-Bessel function of the first kind |
||
real(kind=dp), | intent(inout) | :: | j_ |
j_{l} (x) - Riccati-Bessel funciton of the first kind |
||
real(kind=dp), | intent(inout) | :: | jp_ |
j_{l}' (x) - derivative of the Riccati-Bessel funciton of the first kind |
calculates the Riccati-Bessel function of the second kind and its first derivative. Calls the rcty function from special_functions.f90
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | l_ |
l - order of the Riccati-Bessel function of the second kind |
||
real(kind=dp), | intent(in) | :: | x_ |
x - argument of the Riccati-Bessel function of the second kind |
||
real(kind=dp), | intent(inout) | :: | y_ |
y_{l} (x) - Riccati-Bessel funciton of the second kind |
||
real(kind=dp), | intent(inout) | :: | yp_ |
y_{l}' (x) - derivative of the Riccati-Bessel funciton of the second kind |
determines b, c and d coefficients of the cubic spline function y(x) = y_i + b_i * dx + c_i * dx^2 + d_i * dx^3, where dx = x - x_i, and x_i <= x < x_i+1. The algorithm is based on Gerald, C., and Wheatley, P., "Applied Numerical Analysis", Addison-Wesley, 1994.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | N_ |
number of grid points |
||
real(kind=dp), | intent(in) | :: | x_(N_) |
grid points (ascending order) |
||
real(kind=dp), | intent(in) | :: | y_(N_) |
tabulated values |
||
real(kind=dp), | intent(out) | :: | b_(N_) |
arrays with coefficients of the spline function |
||
real(kind=dp), | intent(out) | :: | c_(N_) |
arrays with coefficients of the spline function |
||
real(kind=dp), | intent(out) | :: | d_(N_) |
arrays with coefficients of the spline function |