The Elk Code
 
Loading...
Searching...
No Matches
sbesseldm.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: sbesseldm
8! !INTERFACE:
9subroutine sbesseldm(m,lmax,x,djl)
10! !INPUT/OUTPUT PARAMETERS:
11! m : order of derivatve (in,integer)
12! lmax : maximum order of Bessel function (in,integer)
13! x : real argument (in,real)
14! djl : array of returned values (out,real(0:lmax))
15! !DESCRIPTION:
16! Computes the $m$th derivative of the spherical Bessel function of the first
17! kind, $j_l(x)$, for argument $x$ and $l=0,1,\ldots,l_{\rm max}$. For
18! $x\ge 1$ this is done by repeatedly using the relations
19! \begin{align*}
20! \frac{d}{dx}j_l(x)&=\frac{l}{x}j_l(x)-j_{l+1}(x) \\
21! j_{l+1}(x)&=\frac{2l+1}{x}j_l(x)-j_{l-1}(x).
22! \end{align*}
23! While for $x<1$ the series expansion of the Bessel function is used
24! $$ \frac{d^m}{dx^m}j_l(x)=\sum_{i=0}^{\infty}
25! \frac{(2i+l)!}{(-2)^ii!(2i+l-m)!(2i+2l+1)!!}x^{2i+l-m}. $$
26! This procedure is numerically stable and accurate to near machine precision
27! for $l\le 30$ and $m\le 6$.
28!
29! !REVISION HISTORY:
30! Created March 2003 (JKD)
31! Modified to return an array of values, October 2004 (JKD)
32!EOP
33!BOC
34implicit none
35! arguments
36integer, intent(in) :: m,lmax
37real(8), intent(in) :: x
38real(8), intent(out) :: djl(0:lmax)
39! local variables
40integer, parameter :: maxm=6,maxns=20
41integer i,j,l,i0
42real(8) t1,sm,x2
43! automatic arrays
44integer a(0:maxm+1),a1(0:maxm+1)
45integer b(0:maxm+1),b1(0:maxm+1)
46real(8) jl(0:lmax+1)
47! external functions
48real(8), external :: factn,factn2,factr
49if (m == 0) then
50 call sbessel(lmax,x,djl)
51 return
52end if
53if ((m < 0).or.(m > maxm)) then
54 write(*,*)
55 write(*,'("Error(sbesseldm): m out of range : ",I8)') m
56 write(*,*)
57 stop
58end if
59if ((lmax < 0).or.(lmax > 30)) then
60 write(*,*)
61 write(*,'("Error(sbesseldm): lmax out of range : ",I8)') lmax
62 write(*,*)
63 stop
64end if
65if ((x < 0.d0).or.(x > 1.d5)) then
66 write(*,*)
67 write(*,'("Error(sbesseldm): x out of range : ",G18.10)') x
68 write(*,*)
69 stop
70end if
71if (x > 1.d0) then
72 call sbessel(lmax+1,x,jl)
73 do l=0,lmax
74 a(1:m+1)=0
75 a(0)=1
76 a1(0:m+1)=0
77 do i=1,m
78 b(0)=0
79 b1(0)=0
80 do j=0,i
81 b(j+1)=a(j)*(l-j)
82 b1(j+1)=-a1(j)*(j+l+2)
83 end do
84 do j=0,i
85 b1(j)=b1(j)-a(j)
86 b(j)=b(j)+a1(j)
87 end do
88 a(0:i+1)=b(0:i+1)
89 a1(0:i+1)=b1(0:i+1)
90 end do
91 t1=1.d0
92 sm=dble(a(0))*jl(l)+dble(a1(0))*jl(l+1)
93 do i=1,m+1
94 t1=t1*x
95 sm=sm+(dble(a(i))*jl(l)+dble(a1(i))*jl(l+1))/t1
96 end do
97 djl(l)=sm
98 end do
99else
100 x2=x**2
101 do l=0,lmax
102 i0=max((m-l+1)/2,0)
103 j=2*i0+l-m
104 if (j == 0) then
105 t1=1.d0
106 else
107 t1=x**j
108 end if
109 t1=factr(j+m,j)*t1/(factn(i0)*factn2(j+l+m+1)*dble((-2)**i0))
110 sm=t1
111 do i=i0+1,maxns
112 j=2*i+l
113 t1=-t1*dble((j-1)*j)*x2/dble((j-l)*(j-m-1)*(j-m)*(j+l+1))
114 if (abs(t1) < 1.d-40) exit
115 sm=sm+t1
116 end do
117 djl(l)=sm
118 end do
119end if
120end subroutine
121!EOC
122
elemental real(8) function factn2(n)
Definition factn2.f90:7
elemental real(8) function factn(n)
Definition factn.f90:7
real(8) function factr(n, d)
Definition factr.f90:10
subroutine sbessel(lmax, x, jl)
Definition sbessel.f90:10
subroutine sbesseldm(m, lmax, x, djl)
Definition sbesseldm.f90:10