src/math.f90: use the intrinsic gamma function

- remove the `gamm` function
- rewrite `fact` in terms of `gamma`
This commit is contained in:
Michele Guerini Rocco 2021-12-20 18:47:28 +01:00
parent 98599b2b7d
commit a30cd09e8e
Signed by: rnhmjoj
GPG Key ID: BFBAF4C975F76450
2 changed files with 10 additions and 36 deletions

View File

@ -904,7 +904,6 @@ end subroutine antihermitian
! !
! !
subroutine ssbi(zz,n,l,fsbi) subroutine ssbi(zz,n,l,fsbi)
use math, only : gamm
implicit none implicit none
! local constants ! local constants
integer, parameter :: lmx=20,nmx=lmx+2 integer, parameter :: lmx=20,nmx=lmx+2
@ -918,7 +917,7 @@ subroutine ssbi(zz,n,l,fsbi)
real(wp_) :: c0,c1,sbi real(wp_) :: c0,c1,sbi
! !
do m=n,l+2 do m=n,l+2
c0=one/gamm(dble(m)+1.5_wp_) c0=one/gamma(dble(m)+1.5_wp_)
sbi=c0 sbi=c0
do k=1,50 do k=1,50
c1=c0*0.25_wp_*zz**2/(dble(m+k)+0.5_wp_)/dble(k) c1=c0*0.25_wp_*zz**2/(dble(m+k)+0.5_wp_)/dble(k)

View File

@ -84,42 +84,17 @@ contains
end if end if
end function catand end function catand
function fact(k) function fact(k)
! Factorial function k!
! Note: the result is real
implicit none implicit none
integer, intent(in) :: k integer, intent(in) :: k
real(wp_) :: fact real(wp_) :: fact
integer :: i
! Factorial function
fact=zero
if(k<0) return
fact=one
if(k==0) return
do i=1,k
fact=fact*i
end do
end function fact
function gamm(xx) fact = gamma(real(1 + k, kind=wp_))
implicit none
real(wp_) :: gamm
real(wp_), intent(in) :: xx
! Returns the value Gamma(xx) for xx > 0.
INTEGER :: j
real(wp_) :: ser,tmp,x,y
real(wp_), parameter :: stp=2.5066282746310005_wp_
real(wp_), dimension(6), parameter :: cof=(/76.18009172947146_wp_, &
-86.50532032941677_wp_,24.01409824083091_wp_,-1.231739572450155_wp_, &
.1208650973866179e-2_wp_,-.5395239384953e-5_wp_/)
x=xx
y=x
tmp=x+5.5_wp_
tmp=(x+0.5_wp_)*log(tmp)-tmp
ser=1.000000000190015_wp_
do j=1,6
y=y+1._wp_
ser=ser+cof(j)/y
end do
gamm=exp(tmp)*(stp*ser/x)
end function gamm
end module math end function fact
end module math