************80
! ITSL0 integrates the Struve function L0(t) from 0 to x.
Discussion:
This procedure evaluates the integral of modified Struve function
L0(t) with respect to t from 0 to x.
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:
31 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, real ( kind = 8 ) X, the upper limit of the integral.
Output, real ( kind = 8 ) TL0, the integral of L0(t) from 0 to x.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=8) | :: | x | ||||
real(kind=8) | :: | tl0 |
subroutine itsl0 ( x, tl0 ) !*****************************************************************************80 ! !! ITSL0 integrates the Struve function L0(t) from 0 to x. ! ! Discussion: ! ! This procedure evaluates the integral of modified Struve function ! L0(t) with respect to t from 0 to x. ! ! 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: ! ! 31 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, real ( kind = 8 ) X, the upper limit of the integral. ! ! Output, real ( kind = 8 ) TL0, the integral of L0(t) from 0 to x. ! implicit none real ( kind = 8 ) a(18) real ( kind = 8 ) a0 real ( kind = 8 ) a1 real ( kind = 8 ) af real ( kind = 8 ) el integer ( kind = 4 ) k real ( kind = 8 ) pi real ( kind = 8 ) r real ( kind = 8 ) rd real ( kind = 8 ) s real ( kind = 8 ) s0 real ( kind = 8 ) ti real ( kind = 8 ) tl0 real ( kind = 8 ) x pi = 3.141592653589793D+00 r = 1.0D+00 if ( x <= 20.0D+00 ) then s = 0.5D+00 do k = 1, 100 if ( k == 1 ) then rd = 0.5D+00 else rd = 1.0D+00 end if r = r * rd * k / ( k + 1.0D+00 ) & * ( x / ( 2.0D+00 * k + 1.0D+00 ) ) ** 2 s = s + r if ( abs ( r / s ) < 1.0D-12 ) then exit end if end do tl0 = 2.0D+00 / pi * x * x * s else s = 1.0D+00 do k = 1, 10 r = r * k / ( k + 1.0D+00 ) & * ( ( 2.0D+00 * k + 1.0D+00 ) / x ) ** 2 s = s + r if ( abs ( r / s ) < 1.0D-12 ) then exit end if end do el = 0.57721566490153D+00 s0 = - s / ( pi * x * x ) + 2.0D+00 / pi & * ( log ( 2.0D+00 * x ) + el ) a0 = 1.0D+00 a1 = 5.0D+00 / 8.0D+00 a(1) = a1 do k = 1, 10 af = ( ( 1.5D+00 * ( k + 0.50D+00 ) & * ( k + 5.0D+00 / 6.0D+00 ) * a1 - 0.5D+00 & * ( k + 0.5D+00 ) ** 2 * ( k -0.5D+00 ) * a0 ) ) & / ( k + 1.0D+00 ) a(k+1) = af a0 = a1 a1 = af end do ti = 1.0D+00 r = 1.0D+00 do k = 1, 11 r = r / x ti = ti + a(k) * r end do tl0 = ti / sqrt ( 2.0D+00 * pi * x ) * exp ( x ) + s0 end if return end subroutine itsl0