853 lines
27 KiB
Fortran
853 lines
27 KiB
Fortran
module conical
|
|
|
|
use const_and_precisions, only : wp_
|
|
implicit none
|
|
|
|
contains
|
|
|
|
function fconic(x,tau,m)
|
|
!
|
|
! this function subprogram computes the conical functions of the
|
|
! first kind P sub(-1/2 + i*tau) (x) for m = 0 and m = 1.
|
|
! Ref. in Kolbig, Comp. Phys. Commun. 23 (1981) 51
|
|
!
|
|
implicit none
|
|
real(wp_), intent(in) :: x, tau
|
|
integer, intent(in) :: m
|
|
real(wp_) :: fconic
|
|
real(wp_) :: t(7),h(9),v(11)
|
|
real(wp_) :: aa,a0,a1,a2,a3,a4,a5,a6,b0,b1,fm,fn,fn1,r1,r2,s,s0,s1
|
|
real(wp_) :: x1,y,y2,y3,z
|
|
integer :: jp,j,n
|
|
real(wp_), parameter :: rpi=1.7724538509055_wp_,pi2=0.63661977236758_wp_
|
|
real(wp_), parameter :: eps=1.0e-14_wp_
|
|
integer, parameter :: nout=2,nmax=200
|
|
!
|
|
complex(wp_) a,b,c,ti,r,rr,q,u,u0,u1,u2,uu
|
|
complex(wp_) v0,v1,v2,vv,w(19)
|
|
!
|
|
logical lm0,lm1,lta
|
|
|
|
fconic=0.0_wp_
|
|
lm0=m == 0
|
|
lm1=m == 1
|
|
if(.not.(lm0 .or. lm1)) then
|
|
write(nout,"(1x,'fconic ... illegal value for m = ',i4)") m
|
|
return
|
|
end if
|
|
fm=m
|
|
fconic=1.0_wp_-fm
|
|
if(x == 1.0_wp_) return
|
|
!
|
|
fconic=0.0_wp_
|
|
if(tau == 0.0_wp_ .and. abs(x-1.0_wp_) > 0.01_wp_) then
|
|
if(x > 1.0_wp_) then
|
|
y=sqrt((x-1.0_wp_)/(x+1.0_wp_))
|
|
z=ellick(y)
|
|
s=sqrt(0.5_wp_*(x+1.0_wp_))
|
|
if(lm0) fconic=pi2*z/s
|
|
if(lm1) fconic=pi2*s*(ellice(y)-z)/sqrt(x**2-1.0_wp_)
|
|
return
|
|
else
|
|
y=sqrt(0.5_wp_*(1.0_wp_-x))
|
|
z=ellick(y)
|
|
if(lm0) fconic=pi2*z
|
|
if(lm1) fconic=pi2*(ellice(y)-0.5_wp_*(1.0_wp_+x)*z)/ &
|
|
sqrt(1.0_wp_-x**2)
|
|
return
|
|
end if
|
|
else
|
|
ti=cmplx(0._wp_,tau,wp_)
|
|
!
|
|
if((-1._wp_ < x .and. x <= 0.0_wp_).or. &
|
|
(0.0_wp_ < x .and. x <= 0.1_wp_ .and.tau<= 17.0_wp_).or. &
|
|
(0.1_wp_ < x .and. x <= 0.2_wp_ .and.tau<= 5.0_wp_)) then
|
|
lta=tau <= 10.0_wp_
|
|
x1=x**2
|
|
a=0.5_wp_*(0.5_wp_-fm-ti)
|
|
b=0.5_wp_*(0.5_wp_-fm+ti)
|
|
c=0.5_wp_
|
|
jp=30
|
|
else if((0.1_wp_ < x .and. x <= 0.2_wp_ .and.tau<= 17.0_wp_) &
|
|
.or.(0.2_wp_ < x .and. x <= 1.5_wp_ .and.tau<= 20.0_wp_)) &
|
|
then
|
|
lta=x > 1.0_wp_ .or. x <= 1.0_wp_ .and. tau <= 5.0_wp_
|
|
x1=(1.0_wp_-x)/2._wp_
|
|
a=0.5_wp_+fm-ti
|
|
b=0.5_wp_+fm+ti
|
|
c=fm+1.0_wp_
|
|
jp=32
|
|
else if(1.5_wp_ < x .and. tau <= max(20.0_wp_,x)) then
|
|
lta=.true.
|
|
x1=1.0_wp_/x**2
|
|
u=exp((-0.5_wp_+ti)*log(2.0_wp_*x)+clogam(1.0_wp_+ti) &
|
|
-clogam(1.5_wp_-fm+ti))
|
|
a=0.5_wp_*(0.5_wp_-fm-ti)
|
|
b=0.5_wp_*(1.5_wp_-fm-ti)
|
|
c=1.0_wp_-ti
|
|
jp=33
|
|
else
|
|
if(x > 1.0_wp_) then
|
|
s=sqrt(x**2-1.0_wp_)
|
|
t(1)=log(x+s)
|
|
h(1)=tau*t(1)
|
|
b0=besj0l(h(1))
|
|
b1=besj1l(h(1))
|
|
z=1.0_wp_
|
|
else
|
|
s=sqrt(1.0_wp_-x**2)
|
|
t(1)=acos(x)
|
|
h(1)=tau*t(1)
|
|
b0=besi0(h(1))
|
|
b1=besi1(h(1))
|
|
z=-1.0_wp_
|
|
end if
|
|
h(1)=t(1)*x/s
|
|
v(1)=tau
|
|
do j = 2,7
|
|
t(j)=t(j-1)*t(1)
|
|
h(j)=h(j-1)*h(1)
|
|
end do
|
|
do j = 2,11
|
|
v(j)=v(j-1)*v(1)
|
|
end do
|
|
!
|
|
if(lm1) then
|
|
aa=-1.0_wp_
|
|
a0=3.0_wp_*(1.0_wp_-h(1))/(8.0_wp_*t(1))
|
|
a1=(-15.0_wp_*h(2)+6.0_wp_*h(1)+9.0_wp_+z*8.0_wp_*t(2))/ &
|
|
(128.0_wp_*t(2))
|
|
a2=3.0_wp_*(-35.0_wp_*h(3)-15.0_wp_*h(2)+15.0_wp_*h(1)+35.0_wp_ &
|
|
+z*t(2)*(32.0_wp_*h(1)+8.0_wp_))/(1024.0_wp_*t(3))
|
|
a3=(-4725.0_wp_*h(4)-6300.0_wp_*h(3)-3150.0_wp_*h(2)+3780.0_wp_*h(1) &
|
|
+10395.0_wp_-1216.0_wp_*t(4)+z*t(2)*(6000.0_wp_*h(2) &
|
|
+5760.0_wp_*h(1)+1680.0_wp_)) /(32768.0_wp_*t(4))
|
|
a4=7.0_wp_*(-10395.0_wp_*h(5)-23625.0_wp_*h(4)-28350.0_wp_*h(3) &
|
|
-14850.0_wp_*h(2)+19305.0_wp_*h(1)+57915.0_wp_ &
|
|
-t(4)*(6336.0_wp_*h(1)+6080.0_wp_)+z*t(2)*(16800.0_wp_*h(3) &
|
|
+30000.0_wp_*h(2)+25920.0_wp_*h(1)+7920.0_wp_))/ &
|
|
(262144.0_wp_*t(5))
|
|
a5=(-2837835.0_wp_*h(6)-9168390.0_wp_*h(5)-16372125.0_wp_*h(4) &
|
|
-18918900*h(3) -10135125.0_wp_*h(2)+13783770.0_wp_*h(1) &
|
|
+43648605.0_wp_-t(4)*(3044160.0_wp_*h(2)+5588352.0_wp_*h(1) &
|
|
+4213440.0_wp_)+z*t(2)*(5556600.0_wp_*h(4)+14817600.0_wp_*h(3) &
|
|
+20790000.0_wp_*h(2)+17297280.0_wp_*h(1)+5405400.0_wp_ &
|
|
+323072.0_wp_*t(4)))/ (4194304.0_wp_*t(6))
|
|
a6=0.0_wp_
|
|
else
|
|
aa=0.0_wp_
|
|
a0=1.0_wp_
|
|
a1=(h(1)-1.0_wp_)/(8.0_wp_*t(1))
|
|
a2=(9.0_wp_*h(2)+6.0_wp_*h(1)-15.0_wp_-z*8.0_wp_*t(2))/ &
|
|
(128.0_wp_*t(2))
|
|
a3=5.0_wp_*(15.0_wp_*h(3)+27.0_wp_*h(2)+21.0_wp_*h(1)-63.0_wp_ &
|
|
-z*t(2)*(16.0_wp_*h(1)+24.0_wp_))/(1024.0_wp_*t(3))
|
|
a4=7.0_wp_*(525.0_wp_*h(4)+1500.0_wp_*h(3)+2430.0_wp_*h(2) &
|
|
+1980.0_wp_*h(1)-6435.0_wp_+192.0_wp_*t(4)-z*t(2)* &
|
|
(720.0_wp_*h(2)+1600.0_wp_*h(1)+2160.0_wp_))/(32768.0_wp_*t(4))
|
|
a5=21.0_wp_*(2835.0_wp_*h(5)+11025.0_wp_*h(4)+24750.0_wp_*h(3) &
|
|
+38610.0_wp_*h(2)+32175.0_wp_*h(1)-109395.0_wp_+t(4) &
|
|
*(1984.0_wp_*h(1)+4032.0_wp_)-z*t(2) &
|
|
*(4800.0_wp_*h(3)+15120.0_wp_*h(2)+26400.0_wp_*h(1)+34320.0_wp_)) &
|
|
/(262144.0_wp_*t(5))
|
|
a6=11.0_wp_*(218295.0_wp_*h(6)+1071630.0_wp_*h(5)+3009825.0_wp_*h(4) &
|
|
+6142500.0_wp_*h(3)+9398025.0_wp_*h(2)+7936110.0_wp_*h(1) &
|
|
-27776385.0_wp_+t(4)*(254016.0_wp_*h(2) &
|
|
+749952.0_wp_*h(1)+1100736.0_wp_)-z*t(2)*(441000.0_wp_*h(4) &
|
|
+1814400.0_wp_*h(3)+4127760.0_wp_*h(2)+6552000.0_wp_*h(1) &
|
|
+8353800.0_wp_+31232.0_wp_*t(4)))/(4194304.0_wp_*t(6))
|
|
end if
|
|
s0=a0+(-4.0_wp_*a3/t(1)+a4)/v(4)+(-192.0_wp_*a5/t(3) &
|
|
+144.0_wp_*a6/t(2))/v(8)+z*(-a2/v(2)+(-24.0_wp_*a4/t(2) &
|
|
+12.0_wp_*a5/t(1)-a6)/v(6)+(-1920.0_wp_*a6/t(4))/v(10))
|
|
s1=a1/v(1)+(8.0_wp_*(a3/t(2)-a4/t(1))+a5)/v(5)+(384.0_wp_*a5/t(4) &
|
|
-768.0_wp_*a6/t(3))/v(9)+z*(aa*v(1)+(2.0_wp_*a2/t(1)-a3)/v(3) &
|
|
+(48.0_wp_*a4/t(3)-72.0_wp_*a5/t(2) &
|
|
+18.0_wp_*a6/t(1))/v(7)+(3840.0_wp_*a6/t(5))/v(11))
|
|
fconic=sqrt(t(1)/s)*(b0*s0+b1*s1)
|
|
return
|
|
end if
|
|
!
|
|
do
|
|
if(lta) then
|
|
y=-x1
|
|
y2=y**2
|
|
y3=y**3
|
|
w(1)=a+1.0_wp_
|
|
w(2)=a+2.0_wp_
|
|
w(3)=b+1.0_wp_
|
|
w(4)=b+2.0_wp_
|
|
w(5)=c+1.0_wp_
|
|
w(6)=c*w(5)
|
|
w(7)=a+b
|
|
w(8)=a*b
|
|
w(9)=(w(8)/c)*y
|
|
w(10)=w(1)*w(3)
|
|
w(11)=w(2)*w(4)
|
|
w(12)=1.0_wp_+(w(11)/(2.0_wp_*w(5)))*y
|
|
w(13)=w(7)-6.0_wp_
|
|
w(14)=w(7)+6.0_wp_
|
|
w(15)=2.0_wp_-w(8)
|
|
w(16)=w(15)-2.0_wp_*w(7)
|
|
!
|
|
v0=1.0_wp_
|
|
v1=1.0_wp_+(w(10)/(2.0_wp_*c))*y
|
|
v2=w(12)+(w(10)*w(11)/(12.0_wp_*w(6)))*y2
|
|
u0=1.0_wp_
|
|
u1=v1-w(9)
|
|
u2=v2-w(9)*w(12)+(w(8)*w(10)/(2.0_wp_*w(6)))*y2
|
|
!
|
|
r=1.0_wp_
|
|
n=2
|
|
do
|
|
n=n+1
|
|
if(n > nmax) then
|
|
write(nout,200) x,tau,m
|
|
return
|
|
end if
|
|
rr=r
|
|
fn=n
|
|
h(1)=fn-1.0_wp_
|
|
h(2)=fn-2.0_wp_
|
|
h(3)=fn-3.0_wp_
|
|
h(4)=2.0_wp_*fn
|
|
h(5)=h(4)-3.0_wp_
|
|
h(6)=2.0_wp_*h(5)
|
|
h(7)=4.0_wp_*(h(4)-1.0_wp_)*h(5)
|
|
h(8)=8.0_wp_*h(5)**2*(h(4)-5.0_wp_)
|
|
h(9)=3.0_wp_*fn**2
|
|
w(1)=a+h(1)
|
|
w(2)=a+h(2)
|
|
w(3)=b+h(1)
|
|
w(4)=b+h(2)
|
|
w(5)=c+h(1)
|
|
w(6)=c+h(2)
|
|
w(7)=c+h(3)
|
|
w(8)=h(2)-a
|
|
w(9)=h(2)-b
|
|
w(10)=h(1)-c
|
|
w(11)=w(1)*w(3)
|
|
w(12)=w(5)*w(6)
|
|
!
|
|
w(17)=1.0_wp_+((h(9)+w(13)*fn+w(16))/(h(6)*w(5)))*y
|
|
w(18)=-((w(11)*w(10)/h(6)+(h(9)-w(14)*fn+w(15))* &
|
|
w(11)*y/h(7))/w(12))*y
|
|
w(19)=(w(2)*w(11)*w(4)*w(8)*w(9)/(h(8)*w(7)*w(12)))*y3
|
|
vv=w(17)*v2+w(18)*v1+w(19)*v0
|
|
uu=w(17)*u2+w(18)*u1+w(19)*u0
|
|
r=uu/vv
|
|
if(abs(r-rr) < eps) exit
|
|
v0=v1
|
|
v1=v2
|
|
v2=vv
|
|
u0=u1
|
|
u1=u2
|
|
u2=uu
|
|
end do
|
|
else
|
|
r=1.0_wp_
|
|
q=1.0_wp_
|
|
do n = 1,nmax
|
|
fn=n
|
|
fn1=fn-1.0_wp_
|
|
rr=r
|
|
q=q*x1*(a+fn1)*(b+fn1)/((c+fn1)*fn)
|
|
r=r+q
|
|
if(abs(r-rr) < eps) exit
|
|
end do
|
|
if (n > nmax) then
|
|
write(nout,200) x,tau,m
|
|
return
|
|
end if
|
|
end if
|
|
if (jp/=30) exit
|
|
r1=real(r)/abs(exp(clogam(a+0.5_wp_)))**2
|
|
a=0.5_wp_*(1.5_wp_-fm-ti)
|
|
b=0.5_wp_*(1.5_wp_-fm+ti)
|
|
c=1.5_wp_
|
|
jp=31
|
|
end do
|
|
if (jp==31) then
|
|
r2=real(r)/abs(exp(clogam(a-0.5_wp_)))**2
|
|
fconic=rpi*(r1-2.0_wp_*x*r2)
|
|
if(lm1) fconic=(2.0_wp_/sqrt(1.0_wp_-x1))*fconic
|
|
return
|
|
else if (jp==32) then
|
|
fconic=real(r)
|
|
if(.not.lm0) then
|
|
fconic=0.5_wp_*(tau**2+0.25_wp_)*sqrt(abs(x**2-1.0_wp_))*fconic
|
|
if(x > 1.0_wp_) fconic=-fconic
|
|
end if
|
|
return
|
|
else if (jp==33) then
|
|
fconic=2.0_wp_*real(u*r*(0.5_wp_-fm+ti)/ti)/rpi
|
|
if(lm1) fconic=fconic/sqrt(1.0_wp_-x1)
|
|
return
|
|
end if
|
|
end if
|
|
!
|
|
200 format(1x,'fconic ... convergence difficulties for c function, x = ', &
|
|
e12.4,5x,'tau = ',e12.4,5x,'m = ',i5)
|
|
!
|
|
end function fconic
|
|
|
|
function clogam(z)
|
|
!
|
|
implicit none
|
|
complex(wp_) :: clogam
|
|
complex(wp_), intent(in) :: z
|
|
complex(wp_) :: v,h,r
|
|
integer :: i,n
|
|
real(wp_) :: x,t,a,c,d,e,f
|
|
integer, parameter :: nout=2
|
|
real(wp_), parameter :: pi=3.1415926535898_wp_
|
|
real(wp_), dimension(10), parameter :: b= &
|
|
(/+8.3333333333333e-2_wp_, -2.7777777777778e-3_wp_, &
|
|
+7.9365079365079e-4_wp_, -5.9523809523810e-4_wp_, &
|
|
+8.4175084175084e-4_wp_, -1.9175269175269e-3_wp_, &
|
|
+6.4102564102564e-3_wp_, -2.9550653594771e-2_wp_, &
|
|
+1.7964437236883e-1_wp_, -1.3924322169059e+0_wp_/)
|
|
!
|
|
x=real(z)
|
|
t=aimag(z)
|
|
if(-abs(x) == aint(x) .and. t == 0.0_wp_) then
|
|
write(nout,'(1x,f20.2)') x
|
|
clogam=(0.0_wp_,0.0_wp_)
|
|
return
|
|
end if
|
|
f=abs(t)
|
|
v=cmplx(x,f,wp_)
|
|
if(x < 0.0_wp_) v=1.0_wp_-v
|
|
h=(0.0_wp_,0.0_wp_)
|
|
c=real(v)
|
|
if(c < 7.0_wp_) then
|
|
n=6-int(c)
|
|
h=v
|
|
d=aimag(v)
|
|
a=atan2(d,c)
|
|
do i = 1,n
|
|
c=c+1.0_wp_
|
|
v=cmplx(c,d,wp_)
|
|
h=h*v
|
|
a=a+atan2(d,c)
|
|
end do
|
|
h=cmplx(0.5_wp_*log(real(h)**2+aimag(h)**2),a,wp_)
|
|
v=v+1.0_wp_
|
|
end if
|
|
r=1.0_wp_/v**2
|
|
clogam=0.91893853320467_wp_+(v-0.5_wp_)*log(v)-v+(b(1)+r*(b(2)+r*(b(3) &
|
|
+r*(b(4)+r*(b(5)+r*(b(6)+r*(b(7)+r*(b(8)+r*(b(9)+r*b(10)))))))))) &
|
|
/v-h
|
|
if(x < 0.0_wp_) then
|
|
!
|
|
a=aint(x)-1.0_wp_
|
|
c=pi*(x-a)
|
|
d=pi*f
|
|
e=exp(-2.0_wp_*d)
|
|
f=sin(c)
|
|
e=d+0.5_wp_*log(e*f**2+0.25_wp_*(1.0_wp_-e)**2)
|
|
f=atan2(cos(c)*tanh(d),f)-a*pi
|
|
clogam=1.1447298858494_wp_-cmplx(e,f,wp_)-clogam
|
|
!
|
|
end if
|
|
if(t < 0.0_wp_) clogam=conjg(clogam)
|
|
!
|
|
end function clogam
|
|
|
|
function ellick(xk)
|
|
implicit none
|
|
real(wp_), intent(in) :: xk
|
|
real(wp_) :: ellick, ellice
|
|
integer :: i
|
|
real(wp_) :: eta,pa,pb,pc,pd
|
|
real(wp_), dimension(10), parameter :: &
|
|
a=(/9.6573590280856e-2_wp_, 3.0885146271305e-2_wp_, &
|
|
1.4938013532687e-2_wp_, 8.7898018745551e-3_wp_, &
|
|
6.1796274460533e-3_wp_, 6.8479092826245e-3_wp_, &
|
|
9.8489293221769e-3_wp_, 8.0030039806500e-3_wp_, &
|
|
2.2966348983970e-3_wp_, 1.3930878570066e-4_wp_/), &
|
|
b=(/1.2499999999991e-1_wp_, 7.0312499739038e-2_wp_, &
|
|
4.8828041906862e-2_wp_, 3.7377739758624e-2_wp_, &
|
|
3.0124849012899e-2_wp_, 2.3931913323111e-2_wp_, &
|
|
1.5530941631977e-2_wp_, 5.9739042991554e-3_wp_, &
|
|
9.2155463496325e-4_wp_, 2.9700280966556e-5_wp_/), &
|
|
c=(/4.4314718056089e-1_wp_, 5.6805194567559e-2_wp_, &
|
|
2.1831811676130e-2_wp_, 1.1569595745295e-2_wp_, &
|
|
7.5950934225594e-3_wp_, 7.8204040609596e-3_wp_, &
|
|
1.0770635039866e-2_wp_, 8.6384421736041e-3_wp_, &
|
|
2.4685033304607e-3_wp_, 1.4946621757181e-4_wp_/), &
|
|
d=(/2.4999999999990e-1_wp_, 9.3749999721203e-2_wp_, &
|
|
5.8593661255531e-2_wp_, 4.2717890547383e-2_wp_, &
|
|
3.3478943665762e-2_wp_, 2.6145014700314e-2_wp_, &
|
|
1.6804023346363e-2_wp_, 6.4321465864383e-3_wp_, &
|
|
9.8983328462254e-4_wp_, 3.1859195655502e-5_wp_/)
|
|
!
|
|
if(abs(xk) >= 1.0_wp_) then
|
|
ellick=0.0_wp_
|
|
return
|
|
end if
|
|
eta=1.0_wp_-xk**2
|
|
pa=a(10)
|
|
do i = 1,9
|
|
pa=pa*eta+a(10-i)
|
|
end do
|
|
pa=pa*eta
|
|
pb=b(10)
|
|
do i = 1,9
|
|
pb=pb*eta+b(10-i)
|
|
end do
|
|
pb=pb*eta
|
|
ellick=1.3862943611199_wp_+pa-log(eta)*(0.5_wp_+pb)
|
|
return
|
|
!
|
|
entry ellice(xk)
|
|
!
|
|
if (abs(xk) >= 1.0_wp_) then
|
|
if (abs(xk) > 1.0_wp_) then
|
|
ellick=0.0_wp_
|
|
else
|
|
ellick=1.0_wp_
|
|
end if
|
|
return
|
|
end if
|
|
eta=1.0_wp_-xk**2
|
|
pc=c(10)
|
|
do i = 1,9
|
|
pc=pc*eta+c(10-i)
|
|
end do
|
|
pc=pc*eta
|
|
pd=d(10)
|
|
do i = 1,9
|
|
pd=pd*eta+d(10-i)
|
|
end do
|
|
pd=pd*eta
|
|
ellick=1.0_wp_+pc-log(eta)*pd
|
|
end function ellick
|
|
|
|
function besjy(x)
|
|
implicit none
|
|
real(wp_), intent(in) :: x
|
|
real(wp_) :: besjy,besj0l,besj1l
|
|
real(wp_) :: besy0,besy1
|
|
logical :: l
|
|
real(wp_) :: v,f,a,b,p,q
|
|
integer, parameter :: nout=2
|
|
!
|
|
entry besj0l(x)
|
|
!
|
|
l=.true.
|
|
v=abs(x)
|
|
if(v >= 8.0_wp_) go to 4
|
|
8 f=0.0625_wp_*x**2-2.0_wp_
|
|
a = - 0.0000000000000008_wp_
|
|
b = f * a + 0.0000000000000413_wp_
|
|
a = f * b - a - 0.0000000000019438_wp_
|
|
b = f * a - b + 0.0000000000784870_wp_
|
|
a = f * b - a - 0.0000000026792535_wp_
|
|
b = f * a - b + 0.0000000760816359_wp_
|
|
a = f * b - a - 0.0000017619469078_wp_
|
|
b = f * a - b + 0.0000324603288210_wp_
|
|
a = f * b - a - 0.0004606261662063_wp_
|
|
b = f * a - b + 0.0048191800694676_wp_
|
|
a = f * b - a - 0.0348937694114089_wp_
|
|
b = f * a - b + 0.1580671023320973_wp_
|
|
a = f * b - a - 0.3700949938726498_wp_
|
|
b = f * a - b + 0.2651786132033368_wp_
|
|
a = f * b - a - 0.0087234423528522_wp_
|
|
a = f * a - b + 0.3154559429497802_wp_
|
|
besjy=0.5_wp_*(a-b)
|
|
if(l) return
|
|
!
|
|
a = + 0.0000000000000016_wp_
|
|
b = f * a - 0.0000000000000875_wp_
|
|
a = f * b - a + 0.0000000000040263_wp_
|
|
b = f * a - b - 0.0000000001583755_wp_
|
|
a = f * b - a + 0.0000000052487948_wp_
|
|
b = f * a - b - 0.0000001440723327_wp_
|
|
a = f * b - a + 0.0000032065325377_wp_
|
|
b = f * a - b - 0.0000563207914106_wp_
|
|
a = f * b - a + 0.0007531135932578_wp_
|
|
b = f * a - b - 0.0072879624795521_wp_
|
|
a = f * b - a + 0.0471966895957634_wp_
|
|
b = f * a - b - 0.1773020127811436_wp_
|
|
a = f * b - a + 0.2615673462550466_wp_
|
|
b = f * a - b + 0.1790343140771827_wp_
|
|
a = f * b - a - 0.2744743055297453_wp_
|
|
a = f * a - b - 0.0662922264065699_wp_
|
|
besjy=0.636619772367581_wp_*log(x)*besjy+0.5_wp_*(a-b)
|
|
return
|
|
!
|
|
4 f=256.0_wp_/x**2-2.0_wp_
|
|
b = + 0.0000000000000007_wp_
|
|
a = f * b - 0.0000000000000051_wp_
|
|
b = f * a - b + 0.0000000000000433_wp_
|
|
a = f * b - a - 0.0000000000004305_wp_
|
|
b = f * a - b + 0.0000000000051683_wp_
|
|
a = f * b - a - 0.0000000000786409_wp_
|
|
b = f * a - b + 0.0000000016306465_wp_
|
|
a = f * b - a - 0.0000000517059454_wp_
|
|
b = f * a - b + 0.0000030751847875_wp_
|
|
a = f * b - a - 0.0005365220468132_wp_
|
|
a = f * a - b + 1.9989206986950373_wp_
|
|
p=a-b
|
|
b = - 0.0000000000000006_wp_
|
|
a = f * b + 0.0000000000000043_wp_
|
|
b = f * a - b - 0.0000000000000334_wp_
|
|
a = f * b - a + 0.0000000000003006_wp_
|
|
b = f * a - b - 0.0000000000032067_wp_
|
|
a = f * b - a + 0.0000000000422012_wp_
|
|
b = f * a - b - 0.0000000007271916_wp_
|
|
a = f * b - a + 0.0000000179724572_wp_
|
|
b = f * a - b - 0.0000007414498411_wp_
|
|
a = f * b - a + 0.0000683851994261_wp_
|
|
a = f * a - b - 0.0311117092106740_wp_
|
|
q=8.0_wp_*(a-b)/v
|
|
f=v-0.785398163397448_wp_
|
|
a=cos(f)
|
|
b=sin(f)
|
|
f=0.398942280401432_wp_/sqrt(v)
|
|
if(l) go to 6
|
|
besjy=f*(q*a+p*b)
|
|
return
|
|
6 besjy=f*(p*a-q*b)
|
|
return
|
|
!
|
|
entry besj1l(x)
|
|
!
|
|
l=.true.
|
|
v=abs(x)
|
|
if(v >= 8.0_wp_) go to 5
|
|
3 f=0.0625_wp_*x**2-2.0_wp_
|
|
b = + 0.0000000000000114_wp_
|
|
a = f * b - 0.0000000000005777_wp_
|
|
b = f * a - b + 0.0000000000252812_wp_
|
|
a = f * b - a - 0.0000000009424213_wp_
|
|
b = f * a - b + 0.0000000294970701_wp_
|
|
a = f * b - a - 0.0000007617587805_wp_
|
|
b = f * a - b + 0.0000158870192399_wp_
|
|
a = f * b - a - 0.0002604443893486_wp_
|
|
b = f * a - b + 0.0032402701826839_wp_
|
|
a = f * b - a - 0.0291755248061542_wp_
|
|
b = f * a - b + 0.1777091172397283_wp_
|
|
a = f * b - a - 0.6614439341345433_wp_
|
|
b = f * a - b + 1.2879940988576776_wp_
|
|
a = f * b - a - 1.1918011605412169_wp_
|
|
a = f * a - b + 1.2967175412105298_wp_
|
|
besjy=0.0625_wp_*(a-b)*x
|
|
if(l) return
|
|
!
|
|
b = - 0.0000000000000244_wp_
|
|
a = f * b + 0.0000000000012114_wp_
|
|
b = f * a - b - 0.0000000000517212_wp_
|
|
a = f * b - a + 0.0000000018754703_wp_
|
|
b = f * a - b - 0.0000000568844004_wp_
|
|
a = f * b - a + 0.0000014166243645_wp_
|
|
b = f * a - b - 0.0000283046401495_wp_
|
|
a = f * b - a + 0.0004404786298671_wp_
|
|
b = f * a - b - 0.0051316411610611_wp_
|
|
a = f * b - a + 0.0423191803533369_wp_
|
|
b = f * a - b - 0.2266249915567549_wp_
|
|
a = f * b - a + 0.6756157807721877_wp_
|
|
b = f * a - b - 0.7672963628866459_wp_
|
|
a = f * b - a - 0.1286973843813500_wp_
|
|
a = f * a - b + 0.0406082117718685_wp_
|
|
besjy=0.636619772367581_wp_*log(x)*besjy-0.636619772367581_wp_/x &
|
|
+0.0625_wp_*(a-b)*x
|
|
return
|
|
!
|
|
5 f=256.0_wp_/x**2-2.0_wp_
|
|
b = - 0.0000000000000007_wp_
|
|
a = f * b + 0.0000000000000055_wp_
|
|
b = f * a - b - 0.0000000000000468_wp_
|
|
a = f * b - a + 0.0000000000004699_wp_
|
|
b = f * a - b - 0.0000000000057049_wp_
|
|
a = f * b - a + 0.0000000000881690_wp_
|
|
b = f * a - b - 0.0000000018718907_wp_
|
|
a = f * b - a + 0.0000000617763396_wp_
|
|
b = f * a - b - 0.0000039872843005_wp_
|
|
a = f * b - a + 0.0008989898330859_wp_
|
|
a = f * a - b + 2.0018060817200274_wp_
|
|
p=a-b
|
|
b = + 0.0000000000000007_wp_
|
|
a = f * b - 0.0000000000000046_wp_
|
|
b = f * a - b + 0.0000000000000360_wp_
|
|
a = f * b - a - 0.0000000000003264_wp_
|
|
b = f * a - b + 0.0000000000035152_wp_
|
|
a = f * b - a - 0.0000000000468636_wp_
|
|
b = f * a - b + 0.0000000008229193_wp_
|
|
a = f * b - a - 0.0000000209597814_wp_
|
|
b = f * a - b + 0.0000009138615258_wp_
|
|
a = f * b - a - 0.0000962772354916_wp_
|
|
a = f * a - b + 0.0935555741390707_wp_
|
|
q=8.0_wp_*(a-b)/v
|
|
f=v-2.356194490192345_wp_
|
|
a=cos(f)
|
|
b=sin(f)
|
|
f=0.398942280401432_wp_/sqrt(v)
|
|
if(l) go to 7
|
|
besjy=f*(q*a+p*b)
|
|
return
|
|
7 besjy=f*(p*a-q*b)
|
|
if(x < 0.0_wp_) besjy=-besjy
|
|
return
|
|
!
|
|
entry besy0(x)
|
|
!
|
|
if(x <= 0.0_wp_) go to 9
|
|
l=.false.
|
|
v=x
|
|
if(v >= 8.0_wp_) go to 4
|
|
go to 8
|
|
entry besy1(x)
|
|
!
|
|
if(x <= 0.0_wp_) go to 9
|
|
l=.false.
|
|
v=x
|
|
if(v >= 8.0_wp_) go to 5
|
|
go to 3
|
|
!
|
|
9 besjy=0.0_wp_
|
|
write(nout,"(1x,'besjy ... non-positive argument x = ',e15.4)") x
|
|
end function besjy
|
|
|
|
function besik(x)
|
|
implicit none
|
|
real(wp_), intent(in) :: x
|
|
real(wp_) :: besik,ebesi0,besi0,ebesi1,besi1,ebesk0,besk0,ebesk1,besk1
|
|
logical :: l,e
|
|
real(wp_) :: v,f,a,b,z
|
|
integer, parameter :: nout=2
|
|
!
|
|
entry ebesi0(x)
|
|
!
|
|
e=.true.
|
|
go to 1
|
|
!
|
|
entry besi0(x)
|
|
!
|
|
e=.false.
|
|
1 l=.true.
|
|
v=abs(x)
|
|
if(v >= 8.0_wp_) go to 4
|
|
8 f=0.0625_wp_*x**2-2.0_wp_
|
|
a = 0.000000000000002_wp_
|
|
b = f * a + 0.000000000000120_wp_
|
|
a = f * b - a + 0.000000000006097_wp_
|
|
b = f * a - b + 0.000000000268828_wp_
|
|
a = f * b - a + 0.000000010169727_wp_
|
|
b = f * a - b + 0.000000326091051_wp_
|
|
a = f * b - a + 0.000008738315497_wp_
|
|
b = f * a - b + 0.000192469359688_wp_
|
|
a = f * b - a + 0.003416331766012_wp_
|
|
b = f * a - b + 0.047718748798174_wp_
|
|
a = f * b - a + 0.509493365439983_wp_
|
|
b = f * a - b + 4.011673760179349_wp_
|
|
a = f * b - a + 22.274819242462231_wp_
|
|
b = f * a - b + 82.489032744024100_wp_
|
|
a = f * b - a + 190.494320172742844_wp_
|
|
a = f * a - b + 255.466879624362167_wp_
|
|
besik=0.5_wp_*(a-b)
|
|
if(l .and. e) besik=exp(-v)*besik
|
|
if(l) return
|
|
!
|
|
a = + 0.000000000000003_wp_
|
|
b = f * a + 0.000000000000159_wp_
|
|
a = f * b - a + 0.000000000007658_wp_
|
|
b = f * a - b + 0.000000000318588_wp_
|
|
a = f * b - a + 0.000000011281211_wp_
|
|
b = f * a - b + 0.000000335195256_wp_
|
|
a = f * b - a + 0.000008216025940_wp_
|
|
b = f * a - b + 0.000162708379043_wp_
|
|
a = f * b - a + 0.002536308188086_wp_
|
|
b = f * a - b + 0.030080722420512_wp_
|
|
a = f * b - a + 0.259084432434900_wp_
|
|
b = f * a - b + 1.511535676029228_wp_
|
|
a = f * b - a + 5.283632866873920_wp_
|
|
b = f * a - b + 8.005368868700334_wp_
|
|
a = f * b - a - 4.563433586448395_wp_
|
|
a = f * a - b - 21.057660177402440_wp_
|
|
besik=-log(0.125_wp_*x)*besik+0.5_wp_*(a-b)
|
|
if(e) besik=exp(x)*besik
|
|
return
|
|
!
|
|
4 f=32.0_wp_/v-2.0_wp_
|
|
b = - 0.000000000000001_wp_
|
|
a = f * b - 0.000000000000001_wp_
|
|
b = f * a - b + 0.000000000000004_wp_
|
|
a = f * b - a + 0.000000000000010_wp_
|
|
b = f * a - b - 0.000000000000024_wp_
|
|
a = f * b - a - 0.000000000000104_wp_
|
|
b = f * a - b + 0.000000000000039_wp_
|
|
a = f * b - a + 0.000000000000966_wp_
|
|
b = f * a - b + 0.000000000001800_wp_
|
|
a = f * b - a - 0.000000000004497_wp_
|
|
b = f * a - b - 0.000000000033127_wp_
|
|
a = f * b - a - 0.000000000078957_wp_
|
|
b = f * a - b + 0.000000000029802_wp_
|
|
a = f * b - a + 0.000000001238425_wp_
|
|
b = f * a - b + 0.000000008513091_wp_
|
|
a = f * b - a + 0.000000056816966_wp_
|
|
b = f * a - b + 0.000000513587727_wp_
|
|
a = f * b - a + 0.000007247591100_wp_
|
|
b = f * a - b + 0.000172700630778_wp_
|
|
a = f * b - a + 0.008445122624921_wp_
|
|
a = f * a - b + 2.016558410917480_wp_
|
|
besik=0.199471140200717_wp_*(a-b)/sqrt(v)
|
|
if(e) return
|
|
besik=exp(v)*besik
|
|
return
|
|
!
|
|
entry ebesi1(x)
|
|
!
|
|
e=.true.
|
|
go to 2
|
|
!
|
|
entry besi1(x)
|
|
!
|
|
e=.false.
|
|
2 l=.true.
|
|
v=abs(x)
|
|
if(v >= 8.0_wp_) go to 3
|
|
7 f=0.0625_wp_*x**2-2.0_wp_
|
|
a = + 0.000000000000001_wp_
|
|
b = f * a + 0.000000000000031_wp_
|
|
a = f * b - a + 0.000000000001679_wp_
|
|
b = f * a - b + 0.000000000079291_wp_
|
|
a = f * b - a + 0.000000003227617_wp_
|
|
b = f * a - b + 0.000000111946285_wp_
|
|
a = f * b - a + 0.000003264138122_wp_
|
|
b = f * a - b + 0.000078756785754_wp_
|
|
a = f * b - a + 0.001543019015627_wp_
|
|
b = f * a - b + 0.023993079147841_wp_
|
|
a = f * b - a + 0.287855511804672_wp_
|
|
b = f * a - b + 2.571459906347755_wp_
|
|
a = f * b - a + 16.334550552522066_wp_
|
|
b = f * a - b + 69.395917633734448_wp_
|
|
a = f * b - a + 181.312616040570265_wp_
|
|
a = f * a - b + 259.890237806477292_wp_
|
|
besik=0.0625_wp_*(a-b)*x
|
|
if(l .and. e) besik=exp(-v)*besik
|
|
if(l) return
|
|
!
|
|
a = + 0.000000000000001_wp_
|
|
b = f * a + 0.000000000000042_wp_
|
|
a = f * b - a + 0.000000000002163_wp_
|
|
b = f * a - b + 0.000000000096660_wp_
|
|
a = f * b - a + 0.000000003696783_wp_
|
|
b = f * a - b + 0.000000119367971_wp_
|
|
a = f * b - a + 0.000003202510692_wp_
|
|
b = f * a - b + 0.000070010627855_wp_
|
|
a = f * b - a + 0.001217056994516_wp_
|
|
b = f * a - b + 0.016300049289816_wp_
|
|
a = f * b - a + 0.161074301656148_wp_
|
|
b = f * a - b + 1.101461993004852_wp_
|
|
a = f * b - a + 4.666387026862842_wp_
|
|
b = f * a - b + 9.361617831395389_wp_
|
|
a = f * b - a - 1.839239224286199_wp_
|
|
a = f * a - b - 26.688095480862668_wp_
|
|
besik=log(0.125_wp_*x)*besik+1.0_wp_/x-0.0625_wp_*(a-b)*x
|
|
if(e) besik=exp(x)*besik
|
|
return
|
|
!
|
|
3 f=32.0_wp_/v-2.0_wp_
|
|
b = + 0.000000000000001_wp_
|
|
a = f * b + 0.000000000000001_wp_
|
|
b = f * a - b - 0.000000000000005_wp_
|
|
a = f * b - a - 0.000000000000010_wp_
|
|
b = f * a - b + 0.000000000000026_wp_
|
|
a = f * b - a + 0.000000000000107_wp_
|
|
b = f * a - b - 0.000000000000053_wp_
|
|
a = f * b - a - 0.000000000001024_wp_
|
|
b = f * a - b - 0.000000000001804_wp_
|
|
a = f * b - a + 0.000000000005103_wp_
|
|
b = f * a - b + 0.000000000035408_wp_
|
|
a = f * b - a + 0.000000000081531_wp_
|
|
b = f * a - b - 0.000000000047563_wp_
|
|
a = f * b - a - 0.000000001401141_wp_
|
|
b = f * a - b - 0.000000009613873_wp_
|
|
a = f * b - a - 0.000000065961142_wp_
|
|
b = f * a - b - 0.000000629724239_wp_
|
|
a = f * b - a - 0.000009732146728_wp_
|
|
b = f * a - b - 0.000277205360764_wp_
|
|
a = f * b - a - 0.024467442963276_wp_
|
|
a = f * a - b + 1.951601204652572_wp_
|
|
besik=0.199471140200717_wp_*(a-b)/sqrt(v)
|
|
if(x < 0.0_wp_) besik=-besik
|
|
if(e) return
|
|
besik=exp(v)*besik
|
|
return
|
|
!
|
|
entry ebesk0 (x)
|
|
!
|
|
e=.true.
|
|
go to 11
|
|
!
|
|
entry besk0(x)
|
|
!
|
|
e=.false.
|
|
11 if(x <= 0.0_wp_) go to 9
|
|
l=.false.
|
|
v=x
|
|
if(x < 5.0_wp_) go to 8
|
|
f=20.0_wp_/x-2.0_wp_
|
|
a = - 0.000000000000002_wp_
|
|
b = f * a + 0.000000000000011_wp_
|
|
a = f * b - a - 0.000000000000079_wp_
|
|
b = f * a - b + 0.000000000000581_wp_
|
|
a = f * b - a - 0.000000000004580_wp_
|
|
b = f * a - b + 0.000000000039044_wp_
|
|
a = f * b - a - 0.000000000364547_wp_
|
|
b = f * a - b + 0.000000003792996_wp_
|
|
a = f * b - a - 0.000000045047338_wp_
|
|
b = f * a - b + 0.000000632575109_wp_
|
|
a = f * b - a - 0.000011106685197_wp_
|
|
b = f * a - b + 0.000269532612763_wp_
|
|
a = f * b - a - 0.011310504646928_wp_
|
|
a = f * a - b + 1.976816348461652_wp_
|
|
besik=0.626657068657750_wp_*(a-b)/sqrt(x)
|
|
if(e) return
|
|
z=besik
|
|
besik=0.0_wp_
|
|
if(x < 180.0_wp_) besik=exp(-x)*z
|
|
return
|
|
!
|
|
entry ebesk1(x)
|
|
!
|
|
e=.true.
|
|
go to 12
|
|
!
|
|
entry besk1(x)
|
|
!
|
|
e=.false.
|
|
12 if(x <= 0.0_wp_) go to 9
|
|
l=.false.
|
|
v=x
|
|
if(x < 5.0_wp_) go to 7
|
|
f=20.0_wp_/x-2.0_wp_
|
|
a = + 0.000000000000002_wp_
|
|
b = f * a - 0.000000000000013_wp_
|
|
a = f * b - a + 0.000000000000089_wp_
|
|
b = f * a - b - 0.000000000000663_wp_
|
|
a = f * b - a + 0.000000000005288_wp_
|
|
b = f * a - b - 0.000000000045757_wp_
|
|
a = f * b - a + 0.000000000435417_wp_
|
|
b = f * a - b - 0.000000004645555_wp_
|
|
a = f * b - a + 0.000000057132218_wp_
|
|
b = f * a - b - 0.000000845172048_wp_
|
|
a = f * b - a + 0.000016185063810_wp_
|
|
b = f * a - b - 0.000468475028167_wp_
|
|
a = f * b - a + 0.035465291243331_wp_
|
|
a = f * a - b + 2.071901717544716_wp_
|
|
besik=0.626657068657750_wp_*(a-b)/sqrt(x)
|
|
if(e) return
|
|
z=besik
|
|
besik=0.0_wp_
|
|
if(x < 180.0_wp_) besik=exp(-x)*z
|
|
return
|
|
9 besik=0.0_wp_
|
|
write(nout,"(1x,'besik ... non-positive argument x = ',e15.4)") x
|
|
end function besik
|
|
!
|
|
! routines for conical function: end
|
|
!
|
|
end module conical |