rctj Subroutine

subroutine rctj(n, x, nm, rj, dj)

************80

! RCTJ computes Riccati-Bessel function of the first kind, and derivatives.

Licensing:

This routine is copyrighted by Shanjie Zhang and Jianming Jin.  However,
they give permission to incorporate this routine into a user program
provided that the copyright is acknowledged.

Modified:

18 July 2012

Author:

Shanjie Zhang, Jianming Jin

Reference:

Shanjie Zhang, Jianming Jin,
Computation of Special Functions,
Wiley, 1996,
ISBN: 0-471-11963-6,
LC: QA351.C45.

Parameters:

Input, integer ( kind = 4 ) N, the order of jn(x).

Input, real ( kind = 8 ) X, the argument.

Output, integer ( kind = 4 ) NM, the highest order computed.

Output, real ( kind = 8 ) RJ(0:N), the values of x jn(x).

Output, real ( kind = 8 ) DJ(0:N), the values of [x jn(x)]'.

Arguments

Type IntentOptional Attributes Name
integer(kind=4) :: n
real(kind=8) :: x
integer(kind=4) :: nm
real(kind=8) :: rj(0:n)
real(kind=8) :: dj(0:n)

Calls

proc~~rctj~2~~CallsGraph proc~rctj~2 rctj msta1 msta1 proc~rctj~2->msta1 msta2 msta2 proc~rctj~2->msta2

Source Code

subroutine rctj ( n, x, nm, rj, dj )

  !*****************************************************************************80
  !
  !! RCTJ computes Riccati-Bessel function of the first kind, and derivatives.
  !
  !  Licensing:
  !
  !    This routine is copyrighted by Shanjie Zhang and Jianming Jin.  However, 
  !    they give permission to incorporate this routine into a user program 
  !    provided that the copyright is acknowledged.
  !
  !  Modified:
  !
  !    18 July 2012
  !
  !  Author:
  !
  !    Shanjie Zhang, Jianming Jin
  !
  !  Reference:
  !
  !    Shanjie Zhang, Jianming Jin,
  !    Computation of Special Functions,
  !    Wiley, 1996,
  !    ISBN: 0-471-11963-6,
  !    LC: QA351.C45.
  !
  !  Parameters:
  !
  !    Input, integer ( kind = 4 ) N, the order of jn(x).
  !
  !    Input, real ( kind = 8 ) X, the argument.
  !
  !    Output, integer ( kind = 4 ) NM, the highest order computed.
  !
  !    Output, real ( kind = 8 ) RJ(0:N), the values of x jn(x).
  !
  !    Output, real ( kind = 8 ) DJ(0:N), the values of [x jn(x)]'.
  !
  implicit none

  integer ( kind = 4 ) n

  real ( kind = 8 ) cs
  real ( kind = 8 ) dj(0:n)
  real ( kind = 8 ) f
  real ( kind = 8 ) f0
  real ( kind = 8 ) f1
  integer ( kind = 4 ) k
  integer ( kind = 4 ) m
    ! integer ( kind = 4 ) msta1
  ! integer ( kind = 4 ) msta2
  integer ( kind = 4 ) nm
  real ( kind = 8 ) rj(0:n)
  real ( kind = 8 ) rj0
  real ( kind = 8 ) rj1
  real ( kind = 8 ) x

  nm = n

  if ( abs ( x ) < 1.0D-100 ) then
     do k = 0, n
        rj(k) = 0.0D+00
        dj(k) = 0.0D+00
     end do
     dj(0) = 1.0D+00
     return
  end if

  rj(0) = sin ( x )
  rj(1) = rj(0) / x - cos ( x )
  rj0 = rj(0)
  rj1 = rj(1)

  if ( 2 <= n ) then

     m = msta1 ( x, 200 )

     if ( m < n ) then
        nm = m
     else
        m = msta2 ( x, n, 15 )
     end if

     f0 = 0.0D+00
     f1 = 1.0D-100
     do k = m, 0, -1
        f = ( 2.0D+00 * k + 3.0D+00 ) * f1 / x - f0
        if ( k <= nm ) then
           rj(k) = f
        end if
        f0 = f1
        f1 = f
     end do

     if ( abs ( rj1 ) < abs ( rj0 ) ) then
        cs = rj0 / f
     else
        cs = rj1 / f0
     end if

     do k = 0, nm
        rj(k) = cs * rj(k)
     end do

  end if

  dj(0) = cos ( x )
  do k = 1, nm
     dj(k) = - k * rj(k) / x + rj(k-1)
  end do

  return
end subroutine rctj