eplf/src/eplf_function.irp.f

353 lines
9.0 KiB
Fortran

BEGIN_PROVIDER [ real, eplf_gamma ]
implicit none
BEGIN_DOC
! Value of the gaussian for the EPLF
END_DOC
include 'constants.F'
real :: eps, N
! N = 0.1
! eps = -real(dlog(tiny(1.d0)))
! 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.e5
END_PROVIDER
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
double precision :: ao_eplf_integral
do i=1,ao_num
do j=i,ao_num
ao_eplf_integral_matrix(j,i) = ao_eplf_integral(j,i,dble(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) ]
implicit none
BEGIN_DOC
! Array of all the <phi_i phi_j | exp(-gamma r^2)> for EPLF
END_DOC
integer :: i, j, k, l
double precision :: t
PROVIDE ao_eplf_integral_matrix
PROVIDE mo_coef mo_coef_transp
do i=1,mo_num
do j=i,mo_num
mo_eplf_integral_matrix(j,i) = 0.d0
enddo
do k=1,ao_num
if (mo_coef(k,i) /= 0.) then
do l=1,ao_num
if (abs(ao_eplf_integral_matrix(l,k))>1.d-16) then
t = mo_coef(k,i)*ao_eplf_integral_matrix(l,k)
do j=i,mo_num
mo_eplf_integral_matrix(j,i) += t*mo_coef_transp(j,l)
enddo
endif
enddo
endif
enddo
enddo
do i=1,mo_num
do j=i+1,mo_num
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
integer :: i, j
eplf_up_up = 0.d0
eplf_up_dn = 0.d0
PROVIDE mo_coef_transp
do i=1,mo_closed_num
do j=1,mo_closed_num
double precision :: temp
temp = mo_value_prod_p(i,i)*mo_eplf_integral_matrix(j,j)
eplf_up_up += temp - mo_value_prod_p(j,i)*mo_eplf_integral_matrix(j,i)
eplf_up_dn += temp
enddo
enddo
eplf_up_up *= 2.d0
eplf_up_dn *= 2.d0
integer :: k,l,m
do m=1,two_e_density_num
i=two_e_density_indice(1,m)
j=two_e_density_indice(2,m)
k=two_e_density_indice(3,m)
l=two_e_density_indice(4,m)
temp = mo_value_prod_p(i,j)*mo_eplf_integral_matrix(k,l)
eplf_up_up += two_e_density_value(1,m)*temp
eplf_up_dn += two_e_density_value(2,m)*temp
enddo
END_PROVIDER
BEGIN_PROVIDER [ real, eplf_value_p ]
implicit none
BEGIN_DOC
! Value of the EPLF at the current point.
END_DOC
double precision :: aa, ab
double precision, parameter :: eps = tiny(1.d0)
aa = eplf_up_up
ab = eplf_up_dn
if ( (aa > 0.d0).and.(ab > 0.d0) ) then
aa = min(1.d0,aa)
ab = min(1.d0,ab)
aa = -dlog(aa)
aa = dsqrt(aa)
ab = -dlog(ab)
ab = dsqrt(ab)
eplf_value_p = (aa-ab)/(aa+ab+eps)
else
eplf_value_p = 0.d0
endif
END_PROVIDER
!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
! 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
! dtemp += &
! (x-xa)**i * (x-xb)**j * exp(-(a*(x-xa)**2+b*(x-xb)**2+gmma*(x-xr)**2))
! x = x+dx
! enddo
! ao_eplf_integral_primitive_oneD_numeric = dtemp*dx
!
!end function
!
!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
! double precision :: ao_eplf_integral_primitive_oneD_numeric
! real :: gmma, center(3), c
!
!
! ao_eplf_integral_numeric = 0.d0
! do q=1,ao_prim_num(j)
! do p=1,ao_prim_num(i)
! c = ao_coef(i,p)*ao_coef(j,q)
! integral = &
! ao_eplf_integral_primitive_oneD_numeric( &
! ao_expo(i,p), &
! nucl_coord(ao_nucl(i),1), &
! ao_power(i,1), &
! ao_expo(j,q), &
! nucl_coord(ao_nucl(j),1), &
! ao_power(j,1), &
! gmma, &
! center(1)) * &
! ao_eplf_integral_primitive_oneD_numeric( &
! ao_expo(i,p), &
! nucl_coord(ao_nucl(i),2), &
! ao_power(i,2), &
! ao_expo(j,q), &
! nucl_coord(ao_nucl(j),2), &
! ao_power(j,2), &
! gmma, &
! center(2)) * &
! ao_eplf_integral_primitive_oneD_numeric( &
! ao_expo(i,p), &
! nucl_coord(ao_nucl(i),3), &
! ao_power(i,3), &
! ao_expo(j,q), &
! nucl_coord(ao_nucl(j),3), &
! ao_power(j,3), &
! gmma, &
! center(3))
! ao_eplf_integral_numeric = ao_eplf_integral_numeric + c*integral
! enddo
! enddo
!
!end function
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 ! Exponents
double precision , intent(in) :: gmma ! eplf_gamma
real, intent(in) :: xa,xb,xr ! Centers
integer, intent(in) :: i,j ! Powers of xa and xb
integer :: ii, jj, kk
ASSERT (a>0.)
ASSERT (b>0.)
ASSERT (i>=0)
ASSERT (j>=0)
double precision :: alpha, beta, delta, c
alpha = a*xa*xa + b*xb*xb + gmma*xr*xr
delta = a*xa + b*xb + gmma*xr
beta = a + b + gmma
c = alpha-delta**2/beta
if ( c > 32.d0 ) then ! Cut-off on exp(-32)
ao_eplf_integral_primitive_oneD = 0.d0
return
endif
double precision :: u,v,w
c = exp(-c)
w = 1.d0/sqrt(beta)
u=delta/beta-xa
v=delta/beta-xb
double precision :: fact2, binom, f, ac
ac = 0.d0
integer :: istart(0:1)
istart = (/ 0, 1 /)
do ii=0,i
do jj=istart(mod(ii,2)),j,2
kk=(ii+jj)/2
f=binom(ii,i)*binom(jj,j)*fact2(ii+jj-1)
ac += u**(i-ii)*v**(j-jj)*w**(ii+jj)*f/dble(2**kk)
enddo
enddo
ao_eplf_integral_primitive_oneD = dsqrt_pi*w*c*ac
end function
double precision function ao_eplf_integral(i,j,gmma,center)
implicit none
integer, intent(in) :: i, j
real, intent(in) :: center(3)
double precision, intent(in) :: gmma
!DEC$ ATTRIBUTES FORCEINLINE
integer :: p,q,k
double precision :: integral
double precision :: ao_eplf_integral_primitive_oneD
double precision :: buffer(100)
ASSERT(i>0)
ASSERT(j>0)
ASSERT(ao_prim_num_max < 100)
ASSERT(i<=ao_num)
ASSERT(j<=ao_num)
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( &
ao_expo(i,p), &
nucl_coord(ao_nucl(i),1), &
ao_power(i,1), &
ao_expo(j,q), &
nucl_coord(ao_nucl(j),1), &
ao_power(j,1), &
gmma, &
center(1))
if (integral /= 0.d0) then
integral *= &
ao_eplf_integral_primitive_oneD( &
ao_expo(i,p), &
nucl_coord(ao_nucl(i),2), &
ao_power(i,2), &
ao_expo(j,q), &
nucl_coord(ao_nucl(j),2), &
ao_power(j,2), &
gmma, &
center(2))
if (integral /= 0.d0) then
integral *= &
ao_eplf_integral_primitive_oneD( &
ao_expo(i,p), &
nucl_coord(ao_nucl(i),3), &
ao_power(i,3), &
ao_expo(j,q), &
nucl_coord(ao_nucl(j),3), &
ao_power(j,3), &
gmma, &
center(3))
buffer(p) += integral*ao_coef(i,p)*ao_coef(j,q)
endif
endif
enddo
enddo
do p=1,ao_prim_num_max
ao_eplf_integral += buffer(p)
enddo
end function
double precision function mo_eplf_integral(i,j)
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
mo_eplf_integral += &
mo_coef(k,i)*mo_coef(l,j)*ao_eplf_integral_matrix(k,l)
enddo
endif
enddo
end function