************80
! IKNA compute Bessel function In(x) and Kn(x), 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:
16 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 In(x) and Kn(x).
Input, real ( kind = 8 ) X, the argument.
Output, integer ( kind = 4 ) NM, the highest order computed.
Output, real ( kind = 8 ) BI(0:N), DI(0:N), BK(0:N), DK(0:N),
the values of In(x), In'(x), Kn(x), Kn'(x).
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=4) | :: | n | ||||
real(kind=8) | :: | x | ||||
integer(kind=4) | :: | nm | ||||
real(kind=8) | :: | bi(0:n) | ||||
real(kind=8) | :: | di(0:n) | ||||
real(kind=8) | :: | bk(0:n) | ||||
real(kind=8) | :: | dk(0:n) |
subroutine ikna ( n, x, nm, bi, di, bk, dk ) !*****************************************************************************80 ! !! IKNA compute Bessel function In(x) and Kn(x), 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: ! ! 16 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 In(x) and Kn(x). ! ! Input, real ( kind = 8 ) X, the argument. ! ! Output, integer ( kind = 4 ) NM, the highest order computed. ! ! Output, real ( kind = 8 ) BI(0:N), DI(0:N), BK(0:N), DK(0:N), ! the values of In(x), In'(x), Kn(x), Kn'(x). ! implicit none integer ( kind = 4 ) n real ( kind = 8 ) bi(0:n) real ( kind = 8 ) bi0 real ( kind = 8 ) bi1 real ( kind = 8 ) bk(0:n) real ( kind = 8 ) bk0 real ( kind = 8 ) bk1 real ( kind = 8 ) di(0:n) real ( kind = 8 ) di0 real ( kind = 8 ) di1 real ( kind = 8 ) dk(0:n) real ( kind = 8 ) dk0 real ( kind = 8 ) dk1 real ( kind = 8 ) f real ( kind = 8 ) f0 real ( kind = 8 ) f1 real ( kind = 8 ) g real ( kind = 8 ) g0 real ( kind = 8 ) g1 real ( kind = 8 ) h real ( kind = 8 ) h0 real ( kind = 8 ) h1 integer ( kind = 4 ) k integer ( kind = 4 ) m ! integer ( kind = 4 ) msta1 ! integer ( kind = 4 ) msta2 integer ( kind = 4 ) nm real ( kind = 8 ) s0 real ( kind = 8 ) x nm = n if ( x <= 1.0D-100 ) then do k = 0, n bi(k) = 0.0D+00 di(k) = 0.0D+00 bk(k) = 1.0D+300 dk(k) = -1.0D+300 end do bi(0) = 1.0D+00 di(1) = 0.5D+00 return end if call ik01a ( x, bi0, di0, bi1, di1, bk0, dk0, bk1, dk1 ) bi(0) = bi0 bi(1) = bi1 bk(0) = bk0 bk(1) = bk1 di(0) = di0 di(1) = di1 dk(0) = dk0 dk(1) = dk1 if ( n <= 1 ) then return end if if ( 40.0D+00 < x .and. n < int ( 0.25D+00 * x ) ) then h0 = bi0 h1 = bi1 do k = 2, n h = -2.0D+00 * ( k - 1.0D+00 ) / x * h1 + h0 bi(k) = h h0 = h1 h1 = h end do else 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 + 1.0D+00 ) * f1 / x + f0 if ( k <= nm ) then bi(k) = f end if f0 = f1 f1 = f end do s0 = bi0 / f do k = 0, nm bi(k) = s0 * bi(k) end do end if g0 = bk0 g1 = bk1 do k = 2, nm g = 2.0D+00 * ( k - 1.0D+00 ) / x * g1 + g0 bk(k) = g g0 = g1 g1 = g end do do k = 2, nm di(k) = bi(k-1) - k / x * bi(k) dk(k) = - bk(k-1) - k / x * bk(k) end do return end subroutine ikna