src/dispersion.f90: compute factorials incrementally

This commit is contained in:
Michele Guerini Rocco 2021-12-20 18:46:44 +01:00
parent 91a2e6cf07
commit 98599b2b7d
Signed by: rnhmjoj
GPG Key ID: BFBAF4C975F76450

View File

@ -262,35 +262,36 @@ subroutine warmdisp(xg,yg,mu,npl,nprf,sox,lrm,err,nprr,npri,fast,imx,ex,ey,ez)
end if end if
! !
end subroutine warmdisp end subroutine warmdisp
!
!
!
subroutine diel_tens_fr(yg,mu,npl,a330,epsl,lrm,fast) subroutine diel_tens_fr(yg,mu,npl,a330,epsl,lrm,fast)
! Fully relativistic case computation of dielectric tensor elements ! Fully relativistic case computation of dielectric tensor elements
! up to third order in Larmor radius for hermitian part ! up to third order in Larmor radius for hermitian part
! !
use math, only : fact
implicit none implicit none
! arguments
! subroutine arguments
integer :: lrm,fast integer :: lrm,fast
real(wp_) :: yg,mu,npl real(wp_) :: yg,mu,npl
complex(wp_) :: a330 complex(wp_) :: a330
complex(wp_), dimension(3,3,lrm) :: epsl complex(wp_), dimension(3,3,lrm) :: epsl
! local variables
! local variables
integer :: i,j,l,is,k,lm integer :: i,j,l,is,k,lm
real(wp_) :: cr,ci real(wp_) :: cr,ci
real(wp_) :: npl2,dnl,w,asl,bsl,cmxw,fal real(wp_) :: npl2,dnl,w,asl,bsl,cmxw,fal
real(wp_), dimension(-lrm:lrm,0:2,0:lrm) :: rr real(wp_), dimension(-lrm:lrm,0:2,0:lrm) :: rr
real(wp_), dimension(lrm,0:2,lrm) :: ri real(wp_), dimension(lrm,0:2,lrm) :: ri
complex(wp_) :: ca11,ca12,ca22,ca13,ca23,ca33,cq0p,cq0m,cq1p,cq1m,cq2p complex(wp_) :: ca11,ca12,ca22,ca13,ca23,ca33,cq0p,cq0m,cq1p,cq1m,cq2p
!
npl2=npl**2 npl2=npl**2
dnl=one-npl2 dnl=one-npl2
!
cmxw=one+15.0_wp_/(8.0_wp_*mu)+105.0_wp_/(128.0_wp_*mu**2) cmxw=one+15.0_wp_/(8.0_wp_*mu)+105.0_wp_/(128.0_wp_*mu**2)
cr=-mu*mu/(sqrt_pi*cmxw) cr=-mu*mu/(sqrt_pi*cmxw)
ci=sqrt(2.0_wp_*pi*mu)*mu**2/cmxw ci=sqrt(2.0_wp_*pi*mu)*mu**2/cmxw
!
do l=1,lrm do l=1,lrm
do j=1,3 do j=1,3
do i=1,3 do i=1,3
@ -298,31 +299,50 @@ subroutine diel_tens_fr(yg,mu,npl,a330,epsl,lrm,fast)
end do end do
end do end do
end do end do
!
if (fast<4) then if (fast<4) then
call hermitian(rr,yg,mu,npl,cr,fast,lrm) call hermitian(rr,yg,mu,npl,cr,fast,lrm)
else else
call hermitian_2(rr,yg,mu,npl,cr,fast,lrm) call hermitian_2(rr,yg,mu,npl,cr,fast,lrm)
end if end if
!
call antihermitian(ri,yg,mu,npl,ci,lrm) call antihermitian(ri,yg,mu,npl,ci,lrm)
! block
! compute all factorials incrementally
integer :: fact_l, fact_2l, fact_l2
integer :: fact_lpis, fact_lmis
fact_l = 1 ! l!
fact_2l = 1 ! (2l)!
fact_l2 = 1 ! (l!)²
do l=1,lrm do l=1,lrm
fact_l = fact_l * l
fact_2l = fact_2l * 2*l * (2*l - 1)
fact_l2 = fact_l2 * l**2
lm=l-1 lm=l-1
fal=-0.25_wp_**l*fact(2*l)/(fact(l)**2*yg**(2*lm)) fal=-0.25_wp_**l * fact_2l/(fact_l2 * yg**(2*lm))
ca11=czero ca11=czero
ca12=czero ca12=czero
ca13=czero ca13=czero
ca22=czero ca22=czero
ca23=czero ca23=czero
ca33=czero ca33=czero
fact_lpis = fact_l ! (l+is)!
fact_lmis = fact_l ! (l-is)!
do is=0,l do is=0,l
k=l-is k=l-is
w=(-one)**k w=(-one)**k
!
asl=w/(fact(is+l)*fact(l-is)) asl=w/(fact_lpis * fact_lmis * one)
bsl=asl*(is*is+dble(2*k*lm*(l+is))/(2.0_wp_*l-one)) bsl=asl*(is*is+dble(2*k*lm*(l+is))/(2.0_wp_*l-one))
!
fact_lpis = fact_lpis * (l + is+1)
fact_lmis = merge(fact_lmis / (l - is), 1, is < l)
if(is.gt.0) then if(is.gt.0) then
cq0p=rr(is,0,l)+rr(-is,0,l)+im*ri(is,0,l) cq0p=rr(is,0,l)+rr(-is,0,l)+im*ri(is,0,l)
cq0m=rr(is,0,l)-rr(-is,0,l)+im*ri(is,0,l) cq0m=rr(is,0,l)-rr(-is,0,l)+im*ri(is,0,l)
@ -336,7 +356,7 @@ subroutine diel_tens_fr(yg,mu,npl,a330,epsl,lrm,fast)
cq1m=rr(is,1,l) cq1m=rr(is,1,l)
cq2p=rr(is,2,l) cq2p=rr(is,2,l)
end if end if
!
ca11=ca11+is**2*asl*cq0p ca11=ca11+is**2*asl*cq0p
ca12=ca12+is*l*asl*cq0m ca12=ca12+is*l*asl*cq0m
ca22=ca22+bsl*cq0p ca22=ca22+bsl*cq0p
@ -351,20 +371,20 @@ subroutine diel_tens_fr(yg,mu,npl,a330,epsl,lrm,fast)
epsl(2,3,l) = - im*ca23*fal epsl(2,3,l) = - im*ca23*fal
epsl(3,3,l) = - ca33*fal epsl(3,3,l) = - ca33*fal
end do end do
! end block
cq2p=rr(0,2,0) cq2p=rr(0,2,0)
a330=cq2p a330=cq2p
!
do l=1,lrm do l=1,lrm
epsl(2,1,l) = - epsl(1,2,l) epsl(2,1,l) = - epsl(1,2,l)
epsl(3,1,l) = epsl(1,3,l) epsl(3,1,l) = epsl(1,3,l)
epsl(3,2,l) = - epsl(2,3,l) epsl(3,2,l) = - epsl(2,3,l)
end do end do
!
end subroutine diel_tens_fr end subroutine diel_tens_fr
!
!
!
subroutine hermitian(rr,yg,mu,npl,cr,fast,lrm) subroutine hermitian(rr,yg,mu,npl,cr,fast,lrm)
use eierf, only : calcei3 use eierf, only : calcei3
implicit none implicit none
@ -780,7 +800,6 @@ end function fhermit
! !
! !
subroutine antihermitian(ri,yg,mu,npl,ci,lrm) subroutine antihermitian(ri,yg,mu,npl,ci,lrm)
use math, only : fact
implicit none implicit none
! local constants ! local constants
integer, parameter :: lmx=20,nmx=lmx+2 integer, parameter :: lmx=20,nmx=lmx+2
@ -856,8 +875,16 @@ subroutine antihermitian(ri,yg,mu,npl,ci,lrm)
else else
ee=exp(-mu*(ygn-one+npl*ub)) ee=exp(-mu*(ygn-one+npl*ub))
call ssbi(aa,n,lrm,fsbi) call ssbi(aa,n,lrm,fsbi)
block
! incrementally compute m!
integer :: fact_m
fact_m = 1
do m=n,lrm do m=n,lrm
cm=sqrt_pi*fact(m)*du**(2*m+1) fact_m = fact_m * m
cm=sqrt_pi*fact_m*du**(2*m+1)
cim=0.5_wp_*ci*dnl**m cim=0.5_wp_*ci*dnl**m
mm=m-n+1 mm=m-n+1
fi0m=cm*fsbi(mm) fi0m=cm*fsbi(mm)
@ -867,6 +894,7 @@ subroutine antihermitian(ri,yg,mu,npl,ci,lrm)
ri(n,1,m)=cim*ee*(du*fi1m+ub*fi0m) ri(n,1,m)=cim*ee*(du*fi1m+ub*fi0m)
ri(n,2,m)=cim*ee*(du*du*fi2m+2.0_wp_*du*ub*fi1m+ub*ub*fi0m) ri(n,2,m)=cim*ee*(du*du*fi2m+2.0_wp_*du*ub*fi1m+ub*ub*fi0m)
end do end do
end block
end if end if
end if end if
end do end do
@ -917,51 +945,70 @@ subroutine expinit
end do end do
! !
end subroutine expinit end subroutine expinit
!
!
!
subroutine diel_tens_wr(yg,mu,npl,a330,epsl,lrm) subroutine diel_tens_wr(yg,mu,npl,a330,epsl,lrm)
! Weakly relativistic dielectric tensor computation of dielectric ! Weakly relativistic dielectric tensor computation of dielectric
! tensor elements (Krivenki and Orefice, JPP 30,125 - 1983) ! tensor elements (Krivenki and Orefice, JPP 30,125 - 1983)
! !
use math, only : fact
implicit none implicit none
! arguments
! subroutine arguments
integer :: lrm integer :: lrm
real(wp_) :: yg,npl,mu real(wp_) :: yg,npl,mu
complex(wp_) :: a330,epsl(3,3,lrm) complex(wp_) :: a330,epsl(3,3,lrm)
! local variables
! local variables
integer :: l,lm,is,k integer :: l,lm,is,k
real(wp_) :: npl2,fcl,w,asl,bsl real(wp_) :: npl2,fcl,w,asl,bsl
complex(wp_) :: ca11,ca12,ca13,ca22,ca23,ca33,cq0p,cq0m,cq1p,cq1m,cq2p complex(wp_) :: ca11,ca12,ca13,ca22,ca23,ca33,cq0p,cq0m,cq1p,cq1m,cq2p
complex(wp_), dimension(0:lrm,0:2) :: cefp,cefm complex(wp_), dimension(0:lrm,0:2) :: cefp,cefm
!
npl2=npl*npl npl2=npl*npl
!
call fsup(cefp,cefm,lrm,yg,npl,mu) call fsup(cefp,cefm,lrm,yg,npl,mu)
!
block
! compute all factorials incrementally
integer :: fact_l, fact_2l_l
integer :: fact_lpis, fact_lmis
fact_l = 1 ! l!
fact_2l_l = 1 ! (2l)!/l!
do l=1,lrm do l=1,lrm
fact_l = fact_l * l
fact_2l_l = fact_2l_l * (4*l - 2) ! see http://oeis.org/A001813
lm=l-1 lm=l-1
fcl=0.5_wp_**l*((one/yg)**2/mu)**lm*fact(2*l)/fact(l) fcl=0.5_wp_**l*((one/yg)**2/mu)**lm*fact_2l_l
ca11=czero ca11=czero
ca12=czero ca12=czero
ca13=czero ca13=czero
ca22=czero ca22=czero
ca23=czero ca23=czero
ca33=czero ca33=czero
fact_lpis = fact_l ! (l+is)!
fact_lmis = fact_l ! (l-is)!
do is=0,l do is=0,l
k=l-is k=l-is
w=(-one)**k w=(-one)**k
!
asl=w/(fact(is+l)*fact(l-is)) asl=w/(fact_lpis * fact_lmis * one)
bsl=asl*(is*is+dble(2*k*lm*(l+is))/(2.0_wp_*l-one)) bsl=asl*(is*is+dble(2*k*lm*(l+is))/(2.0_wp_*l-one))
!
fact_lpis = fact_lpis * (l + is+1)
fact_lmis = merge(fact_lmis / (l - is), 1, is < l)
cq0p=mu*cefp(is,0) cq0p=mu*cefp(is,0)
cq0m=mu*cefm(is,0) cq0m=mu*cefm(is,0)
cq1p=mu*npl*(cefp(is,0)-cefp(is,1)) cq1p=mu*npl*(cefp(is,0)-cefp(is,1))
cq1m=mu*npl*(cefm(is,0)-cefm(is,1)) cq1m=mu*npl*(cefm(is,0)-cefm(is,1))
cq2p=cefp(is,1)+mu*npl2*(cefp(is,2)+cefp(is,0)-2.0_wp_*cefp(is,1)) cq2p=cefp(is,1)+mu*npl2*(cefp(is,2)+cefp(is,0)-2.0_wp_*cefp(is,1))
!
ca11=ca11+is**2*asl*cq0p ca11=ca11+is**2*asl*cq0p
ca12=ca12+is*l*asl*cq0m ca12=ca12+is*l*asl*cq0m
ca22=ca22+bsl*cq0p ca22=ca22+bsl*cq0p
@ -976,20 +1023,20 @@ subroutine diel_tens_wr(yg,mu,npl,a330,epsl,lrm)
epsl(2,3,l) = -im*ca23*fcl epsl(2,3,l) = -im*ca23*fcl
epsl(3,3,l) = -ca33*fcl epsl(3,3,l) = -ca33*fcl
end do end do
! end block
cq2p=cefp(0,1)+mu*npl2*(cefp(0,2)+cefp(0,0)-2.0_wp_*cefp(0,1)) cq2p=cefp(0,1)+mu*npl2*(cefp(0,2)+cefp(0,0)-2.0_wp_*cefp(0,1))
a330=-mu*cq2p a330=-mu*cq2p
!
do l=1,lrm do l=1,lrm
epsl(2,1,l) = - epsl(1,2,l) epsl(2,1,l) = - epsl(1,2,l)
epsl(3,1,l) = epsl(1,3,l) epsl(3,1,l) = epsl(1,3,l)
epsl(3,2,l) = - epsl(2,3,l) epsl(3,2,l) = - epsl(2,3,l)
end do end do
!
end subroutine diel_tens_wr end subroutine diel_tens_wr
!
!
!
subroutine fsup(cefp,cefm,lrm,yg,npl,mu) subroutine fsup(cefp,cefm,lrm,yg,npl,mu)
implicit none implicit none
! local constants ! local constants