riccati_bessel_y Subroutine

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


Calls

proc~~riccati_bessel_y~~CallsGraph proc~riccati_bessel_y riccati_bessel_y proc~rcty rcty proc~riccati_bessel_y->proc~rcty proc~integer_to_character integer_to_character proc~riccati_bessel_y->proc~integer_to_character proc~write_warning write_warning proc~riccati_bessel_y->proc~write_warning proc~write_message write_message proc~write_warning->proc~write_message

Called by

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

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
real(kind=dp), public, dimension(l_+1) :: dy_
integer(kind=int32), public :: highest_order_
real(kind=dp), public, dimension(l_+1) :: ry_

Source Code

      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
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: l_
            !! l - order of the Riccati-Bessel function of the second kind
         real(dp), intent(in) :: x_
            !! x - argument of the Riccati-Bessel function of the second kind
         real(dp), intent(inout) :: y_
            !! y_{l} (x) - Riccati-Bessel funciton of the second kind
         real(dp), intent(inout) :: yp_
            !! y_{l}' (x) - derivative of the Riccati-Bessel funciton
            !! of the second kind
         !---------------------------------------------------------------------!
         integer(int32) :: highest_order_
         real(dp), dimension(l_+1) :: ry_, dy_
         !---------------------------------------------------------------------!
         if(l_ == 0) then
            call rcty(l_+1, x_, highest_order_, ry_, dy_)
         else
            call rcty(l_, x_, highest_order_, ry_, dy_)
         endif
         y_  = ry_(l_ + 1)
         yp_ = dy_(l_ + 1)
         !---------------------------------------------------------------------!
         if (highest_order_ < l_) then
            !------------------------------------------------------------------!
            call write_warning("riccati_bessel_j: maximum order of " //        &
               "Riccati-Bessel function: "// trim(adjustl(integer_to_character(&
               highest_order_))) // "is smaller than requested order l = " //  &
               trim(adjustl(integer_to_character(l_))) )
            !------------------------------------------------------------------!
            y_  = ry_(highest_order_)
            yp_ = dy_(highest_order_)
         else 
            y_  = ry_(l_ + 1)
            yp_ = dy_(l_ + 1)
         endif
         !---------------------------------------------------------------------!
      end subroutine riccati_bessel_y