modified_bessel_k_ratio Subroutine

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


Calls

proc~~modified_bessel_k_ratio~~CallsGraph proc~modified_bessel_k_ratio modified_bessel_k_ratio proc~modified_bessel_temme_algorithm modified_bessel_temme_algorithm proc~modified_bessel_k_ratio->proc~modified_bessel_temme_algorithm proc~write_error write_error proc~modified_bessel_temme_algorithm->proc~write_error proc~rgamma rgamma proc~modified_bessel_temme_algorithm->proc~rgamma proc~write_message write_message proc~write_error->proc~write_message

Called by

proc~~modified_bessel_k_ratio~~CalledByGraph proc~modified_bessel_k_ratio modified_bessel_k_ratio proc~calculate_k_matrix calculate_k_matrix proc~calculate_k_matrix->proc~modified_bessel_k_ratio program~scattering SCATTERING program~scattering->proc~calculate_k_matrix

Contents


Variables

Type Visibility Attributes Name Initial
real(kind=dp), public :: ck_
real(kind=dp), public :: dk_
real(kind=dp), public :: ek_
real(kind=dp), public :: order_

Source Code

      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.
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: l_
            !! l - order of the function (without the 1/2 factor!)
         real(dp), intent(in) :: x_
            !! x - argument of the function
         real(dp), intent(inout) :: ratio_
            !! ratio of the modified Bessel function of the second kind
            !! to its derivative
         !---------------------------------------------------------------------!
         real(dp) :: order_, ck_, dk_, ek_
         !---------------------------------------------------------------------!
         order_ = real(l_, dp) + 0.5_dp
         call modified_bessel_temme_algorithm(order_, x_, ck_, dk_, ek_)
         ratio_ = dk_/ck_
         !---------------------------------------------------------------------!
      end subroutine modified_bessel_k_ratio