math_utilities_mod Module

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


Uses

  • module~~math_utilities_mod~~UsesGraph module~math_utilities_mod math_utilities_mod iso_fortran_env iso_fortran_env module~math_utilities_mod->iso_fortran_env module~special_functions_mod special_functions_mod module~math_utilities_mod->module~special_functions_mod module~utility_functions_mod utility_functions_mod module~math_utilities_mod->module~utility_functions_mod module~utility_functions_mod->iso_fortran_env

Used by

  • module~~math_utilities_mod~~UsedByGraph module~math_utilities_mod math_utilities_mod module~radial_coupling_terms_mod radial_coupling_terms_mod module~radial_coupling_terms_mod->module~math_utilities_mod module~boundary_conditions_mod boundary_conditions_mod module~boundary_conditions_mod->module~math_utilities_mod module~pes_matrix_mod pes_matrix_mod module~pes_matrix_mod->module~math_utilities_mod module~pes_matrix_mod->module~radial_coupling_terms_mod program~scattering SCATTERING program~scattering->module~radial_coupling_terms_mod program~scattering->module~boundary_conditions_mod program~scattering->module~pes_matrix_mod module~propagator_mod propagator_mod program~scattering->module~propagator_mod module~propagator_mod->module~pes_matrix_mod

Contents


Functions

public function is_sum_even(x, y, z) result(sum_even)

checks if the sum of 3 integers is an even integer

Arguments

Type IntentOptional 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

Return Value logical

(out) result: true/false

public function ispline(u_, N_, x_, y_, b_, c_, d_) result(spl_result)

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

Arguments

Type IntentOptional 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

Return Value real(kind=dp)

interpolated value at u_

public function percival_seaton_coefficient(j_, j_prime_, lambda_, omega_) result(percival_seaton_coefficient_)

calculates Percival-Seaton coefficients (body-fixed variant)

Arguments

Type IntentOptional 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}\)

Return Value real(kind=dp)

(out) result: percival seaton coefficient in the body-fixed frame

public function rgamma(x, odd, even) result(rgamma_val)

Calculates 1/Gamma(1-X); modernized version of Molscat's rgamma function; see: https://github.com/molscat/molscat/blob/36fa8f93a92f851e9d84245dd6a972e2910541c5/source_code/rbesjy.f --------------------------------------------------------------------!

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: x
real(kind=dp), intent(out) :: odd
real(kind=dp), intent(out) :: even

Return Value real(kind=dp)

public function triangle_inequality_holds(x, y, z) result(holds)

check if the triangle inequality for 3 variables hols

Arguments

Type IntentOptional 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

Return Value logical

(out) result: true/false

public function zero_projections_3j_condition(x, y, z) result(is_valid)

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

Arguments

Type IntentOptional 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

Return Value logical

(out) result: true/false if conditions are met


Subroutines

public subroutine modified_bessel_k_ratio(l_, x_, ratio_)

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.

Arguments

Type IntentOptional 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

public subroutine modified_bessel_temme_algorithm(v, x, ck, dk, ek)

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

Arguments

Type IntentOptional 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

public subroutine riccati_bessel_j(l_, x_, j_, jp_)

calculates the Riccati-Bessel function of the first kind and its first derivative. Calls the rctj function from special_functions.f90

Arguments

Type IntentOptional 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

public subroutine riccati_bessel_y(l_, x_, y_, yp_)

calculates the Riccati-Bessel function of the second kind and its first derivative. Calls the rcty function from special_functions.f90

Arguments

Type IntentOptional 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

public subroutine spline(N_, x_, y_, b_, c_, d_)

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.

Arguments

Type IntentOptional 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