9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-12-22 03:23:29 +01:00

added scan functional

This commit is contained in:
Emmanuel Giner LCT 2019-06-13 15:58:04 +02:00
parent c568f0b943
commit 244b130a72
3 changed files with 152 additions and 16 deletions

View File

@ -22,6 +22,7 @@
BEGIN_PROVIDER[double precision, mos_grad_in_r_array,(mo_num,n_points_final_grid,3)] BEGIN_PROVIDER[double precision, mos_grad_in_r_array,(mo_num,n_points_final_grid,3)]
&BEGIN_PROVIDER[double precision, mos_grad_in_r_array_tranp,(3,mo_num,n_points_final_grid)]
implicit none implicit none
BEGIN_DOC BEGIN_DOC
! mos_grad_in_r_array(i,j,k) = value of the kth component of the gradient of ith mo on the jth grid point ! mos_grad_in_r_array(i,j,k) = value of the kth component of the gradient of ith mo on the jth grid point
@ -35,6 +36,35 @@
do m=1,3 do m=1,3
call dgemm('N','N',mo_num,n_points_final_grid,ao_num,1.d0,mo_coef_transp,mo_num,aos_grad_in_r_array(1,1,m),ao_num,0.d0,mos_grad_in_r_array(1,1,m),mo_num) call dgemm('N','N',mo_num,n_points_final_grid,ao_num,1.d0,mo_coef_transp,mo_num,aos_grad_in_r_array(1,1,m),ao_num,0.d0,mos_grad_in_r_array(1,1,m),mo_num)
enddo enddo
integer :: i,j
do i = 1, n_points_final_grid
do j = 1, mo_num
do m = 1, 3
mos_grad_in_r_array_tranp(m,j,i) = mos_grad_in_r_array(j,i,m)
enddo
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [double precision, alpha_dens_kin_in_r, (n_points_final_grid)]
&BEGIN_PROVIDER [double precision, beta_dens_kin_in_r, (n_points_final_grid)]
implicit none
integer :: i,m,j
alpha_dens_kin_in_r = 0.d0
beta_dens_kin_in_r = 0.d0
do i = 1, n_points_final_grid
do j = 1, elec_alpha_num
do m = 1, 3
alpha_dens_kin_in_r(i) += 0.5d0 * mos_grad_in_r_array_tranp(m,j,i)**2.d0
enddo
enddo
do j = 1, elec_beta_num
do m = 1, 3
beta_dens_kin_in_r(i) += 0.5d0 * mos_grad_in_r_array_tranp(m,j,i)**2.d0
enddo
enddo
enddo
END_PROVIDER END_PROVIDER
BEGIN_PROVIDER[double precision, mos_lapl_in_r_array,(mo_num,n_points_final_grid,3)] BEGIN_PROVIDER[double precision, mos_lapl_in_r_array,(mo_num,n_points_final_grid,3)]

View File

@ -35,11 +35,8 @@ double precision function ec_lyp_88(rho,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,gr
include 'constants.include.F' include 'constants.include.F'
! Input variables ! Input variables
double precision, intent(in) :: rho,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,grad_rho_2 double precision, intent(in) :: rho,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,grad_rho_2
! Local variables ! Local variables
double precision :: a,b,c,d,c_f,omega,delta double precision :: a,b,c,d,c_f,omega,delta
double precision :: rho_13,rho_inv_13,rho_83,rho_113,rho_inv_113,denom double precision :: rho_13,rho_inv_13,rho_83,rho_113,rho_inv_113,denom
double precision :: thr,huge_num,rho_inv double precision :: thr,huge_num,rho_inv
@ -47,8 +44,6 @@ double precision function ec_lyp_88(rho,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,gr
double precision :: tmp1,tmp2,tmp3,tmp4 double precision :: tmp1,tmp2,tmp3,tmp4
double precision :: big1,big2,big3 double precision :: big1,big2,big3
! Output variables
! Constants of the LYP correlation functional ! Constants of the LYP correlation functional
@ -57,15 +52,25 @@ double precision function ec_lyp_88(rho,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,gr
c = 0.2533d0 c = 0.2533d0
d = 0.349d0 d = 0.349d0
thr = 1d-10 ec_lyp_88 = 0.d0
thr = 1d-15
huge_num = 1.d0/thr huge_num = 1.d0/thr
if(dabs(rho_a).lt.thr)then
return
endif
if(dabs(rho_b).lt.thr)then
return
endif
if(rho.lt.0.d0)then if(rho.lt.0.d0)then
print*,'pb !! rho.lt.0.d0' print*,'pb !! rho.lt.0.d0'
stop stop
endif endif
rho_13 = rho**(1d0/3d0)
rho_113 = rho**(11d0/3d0) rho_13 = rho**(1.d0/3.d0)
rho_113 = rho**(11.d0/3.d0)
if(dabs(rho_13) < thr) then if(dabs(rho_13) < thr) then
rho_inv_13 = huge_num rho_inv_13 = huge_num
@ -76,13 +81,13 @@ double precision function ec_lyp_88(rho,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,gr
if (dabs(rho_113) < thr) then if (dabs(rho_113) < thr) then
rho_inv_113 = huge_num rho_inv_113 = huge_num
else else
rho_inv_113 = 1d0/rho_113 rho_inv_113 = 1.d0/rho_113
endif endif
if (dabs(rho) < thr) then if (dabs(rho) < thr) then
rho_inv = huge_num rho_inv = huge_num
else else
rho_inv = 1d0/rho rho_inv = 1.d0/rho
endif endif
! Useful quantities to predefine ! Useful quantities to predefine
@ -90,21 +95,21 @@ double precision function ec_lyp_88(rho,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,gr
denom = 1d0/(1d0 + d*rho_inv_13) denom = 1d0/(1d0 + d*rho_inv_13)
omega = rho_inv_113*exp(-c*rho_inv_13)*denom omega = rho_inv_113*exp(-c*rho_inv_13)*denom
delta = c*rho_inv_13 + d*rho_inv_13*denom delta = c*rho_inv_13 + d*rho_inv_13*denom
c_f = 0.3d0*(3d0*pi*pi)**(2d0/3d0) c_f = 0.3d0*(3.d0*pi*pi)**(2.d0/3.d0)
rho_2 = rho *rho rho_2 = rho *rho
rho_a_2 = rho_a*rho_a rho_a_2 = rho_a*rho_a
rho_b_2 = rho_b*rho_b rho_b_2 = rho_b*rho_b
cst_2_113 = 2d0**(11d0/3d0) cst_2_113 = 2.d0**(11.d0/3.d0)
cst_8_3 = 8d0/3d0 cst_8_3 = 8.d0/3.d0
! first term in the equation (2) of Preuss CPL, 1989 ! first term in the equation (2) of Preuss CPL, 1989
big1 = 4d0*denom*rho_a*rho_b*rho_inv big1 = 4.d0*denom*rho_a*rho_b*rho_inv
tmp1 = cst_2_113*c_f*(rho_a**cst_8_3 + rho_b**cst_8_3) tmp1 = cst_2_113*c_f*(rho_a**cst_8_3 + rho_b**cst_8_3)
tmp2 = (47d0/18d0 - 7d0/18d0*delta)*grad_rho_2 tmp2 = (47.d0/18.d0 - 7.d0/18.d0*delta)*grad_rho_2
tmp3 = - (5d0/2d0 - 1.d0/18d0*delta)*(grad_rho_a_2 + grad_rho_b_2) tmp3 = - (5d0/2d0 - 1.d0/18d0*delta)*(grad_rho_a_2 + grad_rho_b_2)
tmp4 = - (delta - 11d0)/9d0*(rho_a*rho_inv*grad_rho_a_2 + rho_b*rho_inv*grad_rho_b_2) tmp4 = - (delta - 11d0)/9d0*(rho_a*rho_inv*grad_rho_a_2 + rho_b*rho_inv*grad_rho_b_2)
big2 = rho_a*rho_b*(tmp1 + tmp2 + tmp3 + tmp4) big2 = rho_a*rho_b*(tmp1 + tmp2 + tmp3 + tmp4)
@ -114,7 +119,6 @@ double precision function ec_lyp_88(rho,rho_a,rho_b,grad_rho_a_2,grad_rho_b_2,gr
tmp3 = grad_rho_a_2*(2d0/3d0*rho_2 - rho_b_2) tmp3 = grad_rho_a_2*(2d0/3d0*rho_2 - rho_b_2)
big3 = tmp1 + tmp2 + tmp3 big3 = tmp1 + tmp2 + tmp3
ec_lyp_88 = -a*big1 -a*b*omega*big2 -a*b*omega*big3 ec_lyp_88 = -a*big1 -a*b*omega*big2 -a*b*omega*big3
end end

View File

@ -0,0 +1,102 @@
double precision function ec_scan(rho_a,rho_b,tau,grad_rho_2)
include 'constants.include.F'
implicit none
double precision, intent(in) :: rho_a,rho_b,tau,grad_rho_2
double precision :: cst_13,cst_23,cst_43,cst_53,rho_inv,cst_18,cst_3pi2
double precision :: thr,nup,ndo,xi,s,spin_d,drho,drho2,rho,inv_1alph,e_c_lsda1,h0
double precision :: rs,t_w,t_unif,ds_xi,alpha,fc_alpha,step_f,cst_1alph,beta_inf
double precision :: c_1c,c_2c,d_c,e_c_ldsa1,h1,phi,t,beta_rs,gama,a,w_1,g_at2,phi_3,e_c_1
double precision :: b_1c,b_2c,b_3c,dx_xi,gc_xi,e_c_lsda0,w_0,g_inf,cx_xi,x_inf,f0,e_c_0
thr = 1.d-12
nup = max(rho_a,thr)
ndo = max(rho_b,thr)
rho = nup + ndo
ec_scan = 0.d0
if((rho).lt.thr)return
! constants ...
rho_inv = 1.d0/rho
cst_13 = 1.d0/3.d0
cst_23 = 2.d0 * cst_13
cst_43 = 4.d0 * cst_13
cst_53 = 5.d0 * cst_13
cst_18 = 1.d0/8.d0
cst_3pi2 = 3.d0 * pi*pi
drho2 = max(grad_rho_2,thr)
drho = dsqrt(drho2)
spin_d = max(nup-ndo,thr)
c_1c = 0.64d0
c_2c = 1.5d0
d_c = 0.7d0
b_1c = 0.0285764d0
b_2c = 0.0889d0
b_3c = 0.125541d0
gama = 0.031091d0
! correlation energy lsda1
call ec_only_lda_sr(0.d0,nup,ndo,e_c_lsda1)
xi = spin_d/rho
rs = (cst_43 * pi * rho)**(-cst_13)
s = drho/( 2.d0 * cst_3pi2**(cst_13) * rho**cst_43 )
t_w = drho2 * cst_18
ds_xi = 0.5d0 * ( (1.d0+xi)**cst_53 + (1.d0 - xi)**cst_53)
t_unif = 0.3d0 * (cst_3pi2)**cst_23 * rho**cst_53*ds_xi
t_unif = max(t_unif,thr)
alpha = (tau - t_w)/t_unif
cst_1alph= 1.d0 - alpha
cst_1alph= max(cst_1alph,thr)
inv_1alph= 1.d0/cst_1alph
phi = 0.5d0 * ( (1.d0+xi)**cst_23 + (1.d0 - xi)**cst_23)
phi_3 = phi*phi*phi
t = (cst_3pi2/16.d0)**cst_13 * s / (phi * rs**0.5d0)
w_1 = dexp(-e_c_lsda1/(gama * phi_3))
a = beta_rs(rs) /(gama * w_1)
g_at2 = 1.d0/(1.d0 + 4.d0 * a*t*t)
h1 = gama * phi_3 * dlog(1.d0 + w_1 * (1.d0 - g_at2))
! interpolation function
fc_alpha = dexp(-c_1c * alpha * inv_1alph) * step_f(cst_1alph) - d_c * dexp(c_2c * inv_1alph) * step_f(-cst_1alph)
! first part of the correlation energy
e_c_1 = e_c_lsda1 + h1
dx_xi = 0.5d0 * ( (1.d0+xi)**cst_43 + (1.d0 - xi)**cst_43)
gc_xi = (1.d0 - 2.3631d0 * (dx_xi - 1.d0) ) * (1.d0 - xi**12.d0)
e_c_lsda0= - b_1c / (1.d0 + b_2c * rs**0.5d0 + b_3c * rs)
w_0 = dexp(-e_c_lsda0/b_1c) - 1.d0
beta_inf = 0.066725d0 * 0.1d0 / 0.1778d0
cx_xi = -3.d0/(4.d0*pi) * (9.d0 * pi/4.d0)**cst_13 * dx_xi
if(dabs(xi).lt.thr)then
x_inf = 0.128026d0
else
x_inf = (cst_3pi2/16.d0)**cst_23 * beta_inf * phi/(cx_xi - f0)
endif
f0 = -0.9d0
g_inf = 1.d0/(1.d0 + 4.d0 * x_inf * s*s)**0.25d0
h0 = b_1c * dlog(1.d0 + w_0 * (1.d0 * g_inf))
e_c_0 = (e_c_lsda0 + h0) * gc_xi
ec_scan = e_c_1 + fc_alpha * (e_c_0 - e_c_1)
if(isnan(ec_scan))then
print*,'isnan(ec_scan)'
print*,'e_c_1 + fc_alpha * e_c_0'
print*, e_c_1 , fc_alpha , e_c_0
stop
endif
end
double precision function step_f(x)
implicit none
double precision, intent(in) :: x
if(x.lt.0.d0)then
step_f = 0.d0
else
step_f = 1.d0
endif
end
double precision function beta_rs(rs)
implicit none
double precision, intent(in) ::rs
beta_rs(rs) = 0.066725d0 * (1.d0 + 0.1d0 * rs)/(1.d0 + 0.1778d0 * rs)
end