Improved EPLF function for large gammas

This commit is contained in:
Anthony Scemama 2010-10-27 17:11:38 +02:00
parent b638782709
commit 42d8ba45ec
3 changed files with 61 additions and 139 deletions

View File

@ -15,6 +15,12 @@
! !
!end function !end function
double precision function binom(k,n)
implicit double precision(a-h,o-z)
binom=fact(n)/(fact(k)*fact(n-k))
end
double precision function fact2(n) double precision function fact2(n)
implicit none implicit none
integer :: n integer :: n
@ -80,7 +86,7 @@ double precision function rintgauss(n)
integer :: n integer :: n
rintgauss = sqrt(pi) rintgauss = dsqrt_pi
if ( n == 0 ) then if ( n == 0 ) then
return return
else if ( n == 1 ) then else if ( n == 1 ) then

View File

@ -1 +1,2 @@
double precision, parameter :: pi=3.14159265359d0 double precision, parameter :: pi=3.14159265359d0
double precision, parameter :: dsqrt_pi=1.77245385091d0

View File

@ -5,9 +5,12 @@ BEGIN_PROVIDER [ real, eplf_gamma ]
END_DOC END_DOC
include 'constants.F' include 'constants.F'
real :: eps, N real :: eps, N
N = 0.1 ! N = 0.1
eps = -real(dlog(tiny(1.d0))) ! eps = -real(dlog(tiny(1.d0)))
eplf_gamma = (4./3.*pi*density_value_p/N)**(2./3.) * eps ! eplf_gamma = (4./3.*pi*density_value_p/N)**(2./3.) * eps
N=0.001
eps = 50.
eplf_gamma = eps**2 *0.5d0*(4./9.*pi*density_value_p * N**(-1./3.))**(2./3.)
!eplf_gamma = 1.e10 !eplf_gamma = 1.e10
!eplf_gamma = 1.e5 !eplf_gamma = 1.e5
END_PROVIDER END_PROVIDER
@ -339,155 +342,67 @@ double precision function ao_eplf_integral_primitive_oneD(a,xa,i,b,xb,j,gmma,xr)
double precision , intent(in) :: gmma ! eplf_gamma double precision , intent(in) :: gmma ! eplf_gamma
real, intent(in) :: xa,xb,xr ! Centers real, intent(in) :: xa,xb,xr ! Centers
integer, intent(in) :: i,j ! Powers of xa and xb integer, intent(in) :: i,j ! Powers of xa and xb
integer :: ii, jj, kk, ll integer :: ii, jj, kk
real :: xp1,xp
real :: p1,p
double precision :: S(0:i+1,0:j+1)
double precision :: inv_p(2), di(max(i,j)), dj(j), c
ASSERT (a>0.) ASSERT (a>0.)
ASSERT (b>0.) ASSERT (b>0.)
ASSERT (i>=0) ASSERT (i>=0)
ASSERT (j>=0) ASSERT (j>=0)
! Gaussian product double precision :: alpha, beta, delta, c
real :: t(2), xab(2), ab(2) alpha = a*xa*xa + b*xb*xb + gmma*xr*xr
inv_p(1) = 1.d0/(a+b) delta = a*xa + b*xb + gmma*xr
p1 = a+b beta = a + b + gmma
ab(1) = a*b c = alpha-delta**2/beta
inv_p(2) = 1.d0/(p1+gmma)
t(1) = (a*xa+b*xb)
xab(1) = xa-xb
xp1 = t(1)*inv_p(1)
p = p1+gmma
ab(2) = p1*gmma
t(2) = (p1*xp1+gmma*xr)
xab(2) = xp1-xr
xp = t(2)*inv_p(2)
c = dble(ab(1)*inv_p(1)*xab(1)**2 + &
ab(2)*inv_p(2)*xab(2)**2)
! double precision, save :: c_accu(2)
! c_accu(1) += c
! c_accu(2) += 1.d0
! print *, c_accu(1)/c_accu(2)
if ( c > 32.d0 ) then ! Cut-off on exp(-32) if ( c > 32.d0 ) then ! Cut-off on exp(-32)
ao_eplf_integral_primitive_oneD = 0.d0 ao_eplf_integral_primitive_oneD = 0.d0
return return
endif endif
double precision :: u,v,w
c = exp(-c) c = exp(-c)
w = 1.d0/sqrt(beta)
! Obara-Saika recursion u=delta/beta-xa
S(0,0) = 1.d0 v=delta/beta-xb
do ii=1,max(i,j) double precision :: fact2, binom, f, ac
di(ii) = 0.5d0*inv_p(2)*dble(ii) ac = 0.d0
enddo integer :: istart(0:1)
istart = (/ 0, 1 /)
xab(1) = xp-xa do ii=0,i
xab(2) = xp-xb do jj=istart(mod(ii,2)),j,2
S(1,0) = xab(1) ! * S(0,0) kk=(ii+jj)/2
if (i>1) then f=binom(ii,i)*binom(jj,j)*fact2(ii+jj-1)
do ii=1,i-1 ac += u**(i-ii)*v**(j-jj)*w**(ii+jj)*f/dble(2**kk)
S(ii+1,0) = xab(1) * S(ii,0) + di(ii)*S(ii-1,0) enddo
enddo enddo
endif ao_eplf_integral_primitive_oneD = dsqrt_pi*w*c*ac
S(0,1) = xab(2) ! * S(0,0) end function
if (j>1) then
do jj=1,j-1 double precision function rint(i,xa,a,j,xb,b,xr,gmma)
S(0,jj+1) = xab(2) * S(0,jj) + di(jj)*S(0,jj-1) implicit double precision(a-h,o-z)
enddo include 'constants.F'
endif beta=a+b+gmma
w=1.d0/dsqrt(beta)
do jj=1,j alpha=a*xa**2+b*xb**2+gmma*xr**2
S(1,jj) = xab(1) * S(0,jj) + di(jj) * S(0,jj-1) delta=xa*a+xb*b+xr*gmma
do ii=2,i u=delta/beta-xa
S(ii,jj) = xab(1) * S(ii-1,jj) + di(ii-1) * S(ii-2,jj) + di(jj) * S(ii-1,jj-1) v=delta/beta-xb
enddo ac=0.d0
enddo do n=0,i
do m=0,j
ao_eplf_integral_primitive_oneD = dsqrt(pi*inv_p(2))*c*S(i,j) if(mod(n+m,2).eq.0)then
k=(n+m)/2
end function f=binom(n,i)*binom(m,j)*fact2(n+m-1)
ac=ac+u**(i-n)*v**(j-m)*w**(n+m)*f/2.d0**k
endif
enddo
enddo
rint=ac*dsqrt(pi)*w*dexp(-alpha+delta**2/beta)
end
!double precision function ao_eplf_integral_primitive_oneD(a,xa,i,b,xb,j,gmma,xr)
! implicit none
! include 'constants.F'
!!
! real, intent(in) :: a,b,gmma ! Exponents
! double precision, intent(in) :: gmma
! real, intent(in) :: xa,xb,xr ! Centers
! integer, intent(in) :: i,j ! Powers of xa and xb
! integer :: ii, jj, kk, ll
! real :: xp1,xp
! real :: p1,p
! double precision :: xpa, xpb
! double precision :: inv_p(2),S00, c
! double precision :: ObaraS
!!
! ASSERT (a>0.)
! ASSERT (b>0.)
! ASSERT (i>=0)
! ASSERT (j>=0)
!!
! ! Inlined Gaussian products (same as call gaussian_product)
! real :: t(2), xab(2), ab(2)
! inv_p(1) = 1.d0/(a+b)
! p1 = a+b
! ab(1) = a*b
! inv_p(2) = 1.d0/(p1+gmma)
! t(1) = (a*xa+b*xb)
! xab(1) = xa-xb
! xp1 = t(1)*inv_p(1)
! p = p1+gmma
! ab(2) = p1*gmma
! t(2) = (p1*xp1+gmma*xr)
! xab(2) = xp1-xr
! xp = t(2)*inv_p(2)
! c = exp(- real(ab(1)*inv_p(1)*xab(1)**2 + &
! ab(2)*inv_p(2)*xab(2)**2) )
!!
! xpa = xp-xa
! xpb = xp-xb
! S00 = sqrt(real(pi*inv_p(2)))*c
! ao_eplf_integral_primitive_oneD = ObaraS(i,j,xpa,xpb,inv_p(2),S00)
!!
!end function
!!
!recursive double precision function ObaraS(i,j,xpa,xpb,inv_p,S00) result(res)
! implicit none
! integer, intent(in) :: i, j
! double precision, intent(in) :: xpa, xpb, inv_p
! double precision,intent(in) :: S00
!!
! if (i == 0) then
! if (j == 0) then
! res = S00
! else ! (j>0)
! res = xpb*ObaraS(0,j-1,xpa,xpb,inv_p,S00)
! if (j>1) then
! res += 0.5d0*dble(j-1)*inv_p*ObaraS(0,j-2,xpa,xpb,inv_p,S00)
! endif
! endif ! (i==0).and.(j>0)
! else ! (i>0)
! if (j==0) then
! res = xpa*ObaraS(i-1,0,xpa,xpb,inv_p,S00)
! if (i>1) then
! res += 0.5d0*dble(i-1)*inv_p*ObaraS(i-2,0,xpa,xpb,inv_p,S00)
! endif
! else ! (i>0).and.(j>0)
! res = xpa * ObaraS(i-1,j,xpa,xpb,inv_p,S00)
! if (i>1) then
! res += 0.5d0*dble(i-1)*inv_p*ObaraS(i-2,j,xpa,xpb,inv_p,S00)
! endif
! res += 0.5d0*dble(j)*inv_p*ObaraS(i-1,j-1,xpa,xpb,inv_p,S00)
! endif ! (i>0).and.(j>0)
! endif ! (i>0)
!!
!end function
double precision function ao_eplf_integral(i,j,gmma,center) double precision function ao_eplf_integral(i,j,gmma,center)
implicit none implicit none