10
0
mirror of https://gitlab.com/scemama/eplf synced 2024-06-02 03:15:31 +02:00
eplf/src/eplf_function.irp.f

530 lines
14 KiB
FortranFixed
Raw Normal View History

BEGIN_PROVIDER [ real, eplf_gamma ]
implicit none
BEGIN_DOC
! Value of the gaussian for the EPLF
END_DOC
2009-10-29 18:57:46 +01:00
include 'constants.F'
2010-04-29 18:03:04 +02:00
real :: eps
2009-09-11 17:35:23 +02:00
eps = -real(dlog(tiny(1.d0)))
2010-04-29 18:03:04 +02:00
!real :: N
!N = 0.1
!eplf_gamma = (4./(3.*N)*pi*density_value_p)**(2./3.) * eps
eplf_gamma = density_value_p * eps
2010-05-05 13:08:36 +02:00
!eplf_gamma = 1.e10
!eplf_gamma = 1.e5
END_PROVIDER
2009-05-15 17:18:39 +02:00
BEGIN_PROVIDER [ double precision, ao_eplf_integral_matrix, (ao_num,ao_num) ]
implicit none
BEGIN_DOC
! Array of all the <chi_i chi_j | exp(-gamma r^2)> for EPLF
END_DOC
integer :: i, j
2009-05-15 17:18:39 +02:00
double precision :: ao_eplf_integral
do i=1,ao_num
do j=i,ao_num
2009-05-15 17:18:39 +02:00
ao_eplf_integral_matrix(j,i) = ao_eplf_integral(j,i,eplf_gamma,point)
ao_eplf_integral_matrix(i,j) = ao_eplf_integral_matrix(j,i)
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [ double precision, mo_eplf_integral_matrix, (mo_num,mo_num) ]
2009-05-15 17:18:39 +02:00
implicit none
BEGIN_DOC
! Array of all the <chi_i chi_j | exp(-gamma r^2)> for EPLF
END_DOC
integer :: i, j, k, l
2009-06-15 11:33:16 +02:00
double precision :: t
2009-05-15 17:43:02 +02:00
PROVIDE ao_eplf_integral_matrix
PROVIDE mo_coef
do i=1,mo_num
do j=i,mo_num
2009-10-12 17:37:07 +02:00
mo_eplf_integral_matrix(j,i) = 0.d0
2009-05-15 17:43:02 +02:00
enddo
2009-06-15 11:33:16 +02:00
do k=1,ao_num
2010-04-29 17:15:07 +02:00
if (mo_coef(k,i) /= 0.) then
2009-06-15 11:33:16 +02:00
do l=1,ao_num
t = mo_coef(k,i)*ao_eplf_integral_matrix(l,k)
2010-04-29 18:03:04 +02:00
do j=i,mo_num
mo_eplf_integral_matrix(j,i) += t*mo_coef_transp(j,l)
enddo
2009-06-15 11:33:16 +02:00
enddo
endif
2009-05-15 17:43:02 +02:00
enddo
2009-06-15 11:33:16 +02:00
enddo
do i=1,mo_num
do j=i+1,mo_num
2009-05-15 17:18:39 +02:00
mo_eplf_integral_matrix(i,j) = mo_eplf_integral_matrix(j,i)
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [ double precision, eplf_up_up ]
&BEGIN_PROVIDER [ double precision, eplf_up_dn ]
implicit none
BEGIN_DOC
! Value of the d_upup and d_updn quantities needed for EPLF
END_DOC
2010-04-29 12:59:30 +02:00
integer :: i, j
2009-10-12 17:37:07 +02:00
eplf_up_up = 0.d0
eplf_up_dn = 0.d0
PROVIDE mo_coef_transp
2010-04-29 17:15:07 +02:00
do j=1,mo_closed_num
do i=1,mo_closed_num
eplf_up_up += 2.d0* mo_value_p(i)* ( &
2009-05-15 17:43:02 +02:00
mo_value_p(i)*mo_eplf_integral_matrix(j,j) - &
mo_value_p(j)*mo_eplf_integral_matrix(i,j) )
eplf_up_dn += 2.d0* mo_value_p(i)*mo_value_p(i)* &
mo_eplf_integral_matrix(j,j)
2009-05-15 17:18:39 +02:00
enddo
enddo
2010-04-29 17:15:07 +02:00
integer :: k,l,m,n,p
integer :: ik,il,jk,jl
2010-04-29 17:15:07 +02:00
double precision :: ckl
double precision :: phase
integer :: exc
PROVIDE det
PROVIDE elec_num_2
2010-04-29 17:15:07 +02:00
do k=1,det_num
do l=1,det_num
ckl = det_coef(k)*det_coef(l)
exc = det_exc(k,l,3)
if ( exc < 0 ) then
2010-05-28 22:08:08 +02:00
phase = -1.0d0
exc = -exc
else
2010-05-28 22:08:08 +02:00
phase = 1.0d0
endif
if ( exc == 0 ) then
! Sum over all alpha-alpha and beta-beta interactions
do p=1,2
do m=1,elec_num_2(p)-mo_closed_num
jk = det(m,k,p)
jl = det(m,l,p)
do n=1,elec_num_2(p)-mo_closed_num
ik = det(n,k,p)
il = det(n,l,p)
eplf_up_up += phase*ckl*mo_value_p(ik)* ( &
mo_value_p(il)*mo_eplf_integral_matrix(jk,jl) - &
mo_value_p(jl)*mo_eplf_integral_matrix(jk,il) )
enddo
enddo
2010-04-29 17:15:07 +02:00
enddo
! Sum over all alpha-beta interactions
do m=1,elec_beta_num-mo_closed_num
jk = det(m,k,2)
jl = det(m,l,2)
do n=1,elec_alpha_num-mo_closed_num
ik = det(n,k,1)
il = det(n,l,1)
eplf_up_dn += phase*ckl * ( mo_value_p(ik)*mo_value_p(il) * mo_eplf_integral_matrix(jk,jl) &
+ mo_value_p(jk)*mo_value_p(jl) * mo_eplf_integral_matrix(ik,il) )
enddo
enddo
else if ( exc == 1 ) then
do p=1,2
if ( det_exc(k,l,p) == 1 ) then
! Sum over only the sigma-sigma interactions involving the excitation
call get_single_excitation(k,l,ik,il,p)
do m=1,elec_num_2(p)-mo_closed_num
jk = det(m,k,p)
jl = det(m,l,p)
eplf_up_up += phase*ckl*mo_value_p(ik)* ( &
mo_value_p(il)*mo_eplf_integral_matrix(jk,jl) - &
mo_value_p(jl)*mo_eplf_integral_matrix(jk,il) )
enddo
! Sum over only the sigma-(sigma_bar) interactions involving the excitation
integer :: p2
p2 = 1+mod(p,2)
do m=1,elec_num_2(p2)-mo_closed_num
jk = det(m,k,p2)
jl = det(m,l,p2)
eplf_up_dn += phase*ckl * ( mo_value_p(ik)*mo_value_p(il) * mo_eplf_integral_matrix(jk,jl) &
+ mo_value_p(jk)*mo_value_p(jl) * mo_eplf_integral_matrix(ik,il) )
enddo
endif
enddo
else if (exc == 2) then
if ( ( det_exc(k,l,1) == 2 ).or.( det_exc(k,l,2) == 2 ) ) then
! Consider only the double excitations of same-spin electrons
if ( det_exc(k,l,1) == 2 ) then
call get_double_excitation(k,l,ik,jk,il,jl,1)
else if ( det_exc(k,l,2) == 2 ) then
call get_double_excitation(k,l,ik,jk,il,jl,2)
endif
eplf_up_up += phase*ckl*mo_value_p(ik)* ( &
mo_value_p(il)*mo_eplf_integral_matrix(jk,jl) - &
mo_value_p(jl)*mo_eplf_integral_matrix(jk,il) )
else if ( det_exc(k,l,1) == 1 ) then
! Consider only the double excitations of opposite-spin electrons
call get_single_excitation(k,l,ik,jk,1)
call get_single_excitation(k,l,il,jl,2)
eplf_up_dn += phase*ckl * ( mo_value_p(ik)*mo_value_p(il) * mo_eplf_integral_matrix(jk,jl) &
+ mo_value_p(jk)*mo_value_p(jl) * mo_eplf_integral_matrix(ik,il) )
endif
endif
enddo
enddo
2010-04-29 17:15:07 +02:00
END_PROVIDER
2009-11-06 00:27:24 +01:00
BEGIN_PROVIDER [ real, eplf_value_p ]
implicit none
BEGIN_DOC
! Value of the EPLF at the current point.
END_DOC
double precision :: aa, ab
2009-09-11 17:35:23 +02:00
double precision, parameter :: eps = tiny(1.d0)
aa = eplf_up_up
ab = eplf_up_dn
2009-09-11 17:35:23 +02:00
if ( (aa > 0.d0).and.(ab > 0.d0) ) then
aa = min(1.d0,aa)
ab = min(1.d0,ab)
aa = -(dlog(aa)/eplf_gamma)
ab = -(dlog(ab)/eplf_gamma)
aa = dsqrt(aa)
ab = dsqrt(ab)
2009-11-06 00:27:24 +01:00
eplf_value_p = (aa-ab)/(aa+ab+eps)
2009-09-11 17:35:23 +02:00
else
2009-11-06 00:27:24 +01:00
eplf_value_p = 0.d0
2009-09-11 17:35:23 +02:00
endif
END_PROVIDER
2009-05-15 17:18:39 +02:00
double precision function ao_eplf_integral_primitive_oneD_numeric(a,xa,i,b,xb,j,gmma,xr)
implicit none
include 'constants.F'
real, intent(in) :: a,b,gmma ! Exponents
real, intent(in) :: xa,xb,xr ! Centers
integer, intent(in) :: i,j ! Powers of xa and xb
2009-06-15 11:33:16 +02:00
integer,parameter :: Npoints=10000
real :: x, xmin, xmax, dx
ASSERT (a>0.)
ASSERT (b>0.)
ASSERT (i>=0)
ASSERT (j>=0)
xmin = min(xa,xb)
xmax = max(xa,xb)
xmin = min(xmin,xr) - 10.
xmax = max(xmax,xr) + 10.
dx = (xmax-xmin)/real(Npoints)
real :: dtemp
dtemp = 0.
x = xmin
integer :: k
do k=1,Npoints
2009-12-21 23:39:00 +01:00
dtemp += &
(x-xa)**i * (x-xb)**j * exp(-(a*(x-xa)**2+b*(x-xb)**2+gmma*(x-xr)**2))
x = x+dx
enddo
2009-05-15 17:18:39 +02:00
ao_eplf_integral_primitive_oneD_numeric = dtemp*dx
end function
2009-05-15 17:18:39 +02:00
double precision function ao_eplf_integral_numeric(i,j,gmma,center)
implicit none
integer, intent(in) :: i, j
integer :: p,q,k
double precision :: integral
2009-05-15 17:18:39 +02:00
double precision :: ao_eplf_integral_primitive_oneD_numeric
2009-06-15 11:33:16 +02:00
real :: gmma, center(3), c
2009-10-12 17:37:07 +02:00
ao_eplf_integral_numeric = 0.d0
do q=1,ao_prim_num(j)
do p=1,ao_prim_num(i)
2009-10-12 17:37:07 +02:00
c = ao_coef(i,p)*ao_coef(j,q)
integral = &
2009-05-15 17:18:39 +02:00
ao_eplf_integral_primitive_oneD_numeric( &
2009-10-12 17:37:07 +02:00
ao_expo(i,p), &
nucl_coord(ao_nucl(i),1), &
ao_power(i,1), &
2009-10-12 17:37:07 +02:00
ao_expo(j,q), &
nucl_coord(ao_nucl(j),1), &
ao_power(j,1), &
gmma, &
center(1)) * &
2009-05-15 17:18:39 +02:00
ao_eplf_integral_primitive_oneD_numeric( &
2009-10-12 17:37:07 +02:00
ao_expo(i,p), &
nucl_coord(ao_nucl(i),2), &
ao_power(i,2), &
2009-10-12 17:37:07 +02:00
ao_expo(j,q), &
nucl_coord(ao_nucl(j),2), &
ao_power(j,2), &
gmma, &
center(2)) * &
2009-05-15 17:18:39 +02:00
ao_eplf_integral_primitive_oneD_numeric( &
2009-10-12 17:37:07 +02:00
ao_expo(i,p), &
nucl_coord(ao_nucl(i),3), &
ao_power(i,3), &
2009-10-12 17:37:07 +02:00
ao_expo(j,q), &
nucl_coord(ao_nucl(j),3), &
ao_power(j,3), &
gmma, &
center(3))
2009-06-15 11:33:16 +02:00
ao_eplf_integral_numeric = ao_eplf_integral_numeric + c*integral
enddo
enddo
end function
2009-12-22 00:35:27 +01:00
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
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 :: S(0:i+1,0:j+1)
double precision :: inv_p(2), di(max(i,j)), dj(j), c
ASSERT (a>0.)
ASSERT (b>0.)
ASSERT (i>=0)
ASSERT (j>=0)
! Gaussian product
! Inlined Gaussian products (same as call gaussian_product)
2009-12-22 00:07:34 +01:00
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) )
2009-12-22 00:35:27 +01:00
S(0,0) = dsqrt(pi*inv_p(2))*c
! Obara-Saika recursion
do ii=1,max(i,j)
di(ii) = 0.5d0*inv_p(2)*dble(ii)
enddo
S(1,0) = (xp-xa) * S(0,0)
if (i>1) then
do ii=1,i-1
S(ii+1,0) = (xp-xa) * S(ii,0) + di(ii)*S(ii-1,0)
enddo
endif
S(0,1) = (xp-xb) * S(0,0)
if (j>1) then
do jj=1,j-1
S(0,jj+1) = (xp-xb) * S(0,jj) + di(jj)*S(0,jj-1)
enddo
endif
do jj=1,j
S(1,jj) = (xp-xa) * S(0,jj) + di(jj) * S(0,jj-1)
do ii=2,i
S(ii,jj) = (xp-xa) * S(ii-1,jj) + di(ii-1) * S(ii-2,jj) + di(jj) * S(ii-1,jj-1)
enddo
enddo
ao_eplf_integral_primitive_oneD = S(i,j)
end function
2009-12-21 23:39:00 +01:00
2009-12-22 00:35:27 +01:00
!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
! 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
2009-05-15 17:18:39 +02:00
double precision function ao_eplf_integral(i,j,gmma,center)
implicit none
integer, intent(in) :: i, j
integer :: p,q,k
double precision :: integral
2009-12-22 00:35:27 +01:00
!DEC$ ATTRIBUTES FORCEINLINE
2009-05-15 17:18:39 +02:00
double precision :: ao_eplf_integral_primitive_oneD
real :: gmma, center(3)
2009-12-22 00:35:27 +01:00
double precision :: buffer(100)
ASSERT(i>0)
ASSERT(j>0)
ASSERT(i<=ao_num)
ASSERT(j<=ao_num)
2009-12-22 00:35:27 +01:00
ao_eplf_integral = 0.d0
do p=1,ao_prim_num_max
buffer(p) = 0.d0
enddo
do q=1,ao_prim_num(j)
do p=1,ao_prim_num(i)
integral = &
ao_eplf_integral_primitive_oneD( &
2009-10-12 17:37:07 +02:00
ao_expo(i,p), &
nucl_coord(ao_nucl(i),1), &
ao_power(i,1), &
2009-10-12 17:37:07 +02:00
ao_expo(j,q), &
nucl_coord(ao_nucl(j),1), &
ao_power(j,1), &
gmma, &
center(1)) * &
ao_eplf_integral_primitive_oneD( &
2009-10-12 17:37:07 +02:00
ao_expo(i,p), &
nucl_coord(ao_nucl(i),2), &
ao_power(i,2), &
2009-10-12 17:37:07 +02:00
ao_expo(j,q), &
nucl_coord(ao_nucl(j),2), &
ao_power(j,2), &
gmma, &
center(2)) * &
ao_eplf_integral_primitive_oneD( &
2009-10-12 17:37:07 +02:00
ao_expo(i,p), &
nucl_coord(ao_nucl(i),3), &
ao_power(i,3), &
2009-10-12 17:37:07 +02:00
ao_expo(j,q), &
nucl_coord(ao_nucl(j),3), &
ao_power(j,3), &
gmma, &
center(3))
2009-12-22 00:35:27 +01:00
! ao_eplf_integral += integral*ao_coef(i,p)*ao_coef(j,q)
buffer(p) += integral*ao_coef(i,p)*ao_coef(j,q)
enddo
enddo
2009-12-22 00:35:27 +01:00
do p=1,ao_prim_num_max
ao_eplf_integral += buffer(p)
enddo
end function
2009-10-12 17:37:07 +02:00
double precision function mo_eplf_integral(i,j)
2009-05-20 17:57:35 +02:00
implicit none
integer :: i, j, k, l
PROVIDE ao_eplf_integral_matrix
PROVIDE mo_coef
mo_eplf_integral = 0.d0
do k=1,ao_num
if (mo_coef(k,i) /= 0.) then
do l=1,ao_num
2009-12-21 23:39:00 +01:00
mo_eplf_integral += &
2009-05-20 17:57:35 +02:00
mo_coef(k,i)*mo_coef(l,j)*ao_eplf_integral_matrix(k,l)
enddo
endif
enddo
2009-10-12 17:37:07 +02:00
end function
2009-05-20 17:57:35 +02:00