TABLE OF CONTENTS
ABINIT/xctfw [ Functions ]
NAME
xctfw
FUNCTION
Add gradient part of the Thomas-Fermi-Weizsacker functional Perrot F.,Phys. Rev. A20,586-594 (1979)
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (JFD,LK) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors .
INPUTS
ndvxcdgr= size of dvxcdgr(npts,ndvxcdgr) ngr2= size of grho2_updn(npts,ngr2) npts= number of points to be computed nspden=number if spin density component (necessarily 1 here) grho2_updn(npts,ngr2)=square of the gradient of the spin-up, and, if nspden==2, spin-down, and total density (Hartree/Bohr**2), only used if gradient corrected functional (option=2,-2,-4 and 4 or beyond) rho_updn(npts,nspden)=spin-up and spin-down density (Hartree/bohr**3) temp= electronic temperature usefxc=1 if free energy fxc is used
SIDE EFFECTS
The following arrays are modified (gradient correction added): dvxcdgr(npts,3)=partial derivative of the XC energy divided by the norm of the gradient exci(npts)=exchange-correlation energy density fxci(npts)=free energy energy density vxci(npts,nspden)=exchange-correlation potential
PARENTS
rhotoxc
CHILDREN
invcb
SOURCE
43 !!$#if defined HAVE_CONFIG_H 44 !!$#include "config.h" 45 !!$#endif 46 47 #if defined HAVE_CONFIG_H 48 #include "config.h" 49 #endif 50 51 #include "abi_common.h" 52 53 subroutine xctfw(temp,exci,fxci,usefxc,rho_updn,vxci,npts,nspden,dvxcdgr,ndvxcdgr,grho2_updn,ngr2) 54 55 use defs_basis 56 use m_profiling_abi 57 use m_errors 58 59 !This section has been created automatically by the script Abilint (TD). 60 !Do not modify the following lines by hand. 61 #undef ABI_FUNC 62 #define ABI_FUNC 'xctfw' 63 use interfaces_41_xc_lowlevel, except_this_one => xctfw 64 !End of the abilint section 65 66 implicit none 67 68 !Arguments ------------------------------------ 69 !scalars 70 integer,intent(in) :: ndvxcdgr,ngr2,npts,nspden,usefxc 71 real(dp),intent(in) :: temp 72 !arrays 73 real(dp),intent(in) :: grho2_updn(npts,ngr2),rho_updn(npts,nspden) 74 real(dp),intent(inout) :: dvxcdgr(npts,ndvxcdgr),fxci(npts*usefxc),exci(npts),vxci(npts,nspden) 75 76 !Local variables------------------------------- 77 !scalars 78 integer :: iperrot,ipts 79 logical :: has_fxc,has_dvxcdgr 80 real(dp) :: etfw,rho,rho_inv,rhomot,yperrot0,vtfw 81 real(dp) :: yperrot,uperrot,dyperrotdn,duperrotdyperrot 82 real(dp) :: hperrot,dhperrotdyperrot,dhperrotdn,dhperrotduperrot 83 !arrays 84 real(dp) :: wpy(0:7), wpu(0:7) 85 real(dp),allocatable :: rho_updnm1_3(:,:) 86 87 ! ************************************************************************* 88 89 !DBG_ENTER('COLL') 90 91 has_fxc=(usefxc/=0) 92 has_dvxcdgr=(ndvxcdgr/=0) 93 94 yperrot0=1.666081101_dp 95 96 wpy(0)=0.5_dp; wpy(1)=-0.1999176316_dp 97 wpy(2)=0.09765615709_dp; wpy(3)=-0.06237609924_dp 98 wpy(4)=0.05801466322_dp; wpy(5)=-0.04449287774_dp 99 wpy(6)=0.01903211697_dp; wpy(7)=-0.003284096926_dp 100 101 wpu(0)=one/6._dp; wpu(1)=0.311590799_dp 102 wpu(2)=3.295662439_dp; wpu(3)=-29.22038326_dp 103 wpu(4)=116.1084531_dp; wpu(5)=-250.4543147_dp 104 wpu(6)=281.433688_dp; wpu(7)=-128.8784806_dp 105 106 ABI_ALLOCATE(rho_updnm1_3,(npts,2)) 107 108 call invcb(rho_updn(:,1),rho_updnm1_3(:,1),npts) 109 110 do ipts=1,npts 111 112 rho =rho_updn(ipts,1) 113 rhomot=rho_updnm1_3(ipts,1) 114 rho_inv=rhomot*rhomot*rhomot 115 116 yperrot=pi*pi/sqrt2/temp**1.5*two*rho 117 uperrot=yperrot**(2./3.) 118 119 dyperrotdn=pi*pi/sqrt2/temp**1.5*2.0_dp 120 121 hperrot=zero 122 dhperrotdyperrot=zero 123 dhperrotduperrot=zero 124 if(yperrot<=yperrot0)then 125 do iperrot=0,7 126 hperrot=hperrot+wpy(iperrot)*yperrot**iperrot 127 dhperrotdyperrot=dhperrotdyperrot+iperrot*wpy(iperrot)*yperrot**(iperrot-1) 128 end do 129 hperrot=one/12.0_dp*hperrot 130 dhperrotdyperrot=one/12.0_dp*dhperrotdyperrot 131 dhperrotdn=dhperrotdyperrot*dyperrotdn 132 else 133 do iperrot=0,7 134 hperrot=hperrot+wpu(iperrot)/uperrot**(2*iperrot) 135 dhperrotduperrot=dhperrotduperrot-2.*iperrot*wpu(iperrot)/uperrot**(2*iperrot+1) 136 end do 137 hperrot=one/12.0_dp*hperrot 138 dhperrotduperrot=one/12.0_dp*dhperrotduperrot 139 duperrotdyperrot=two/3._dp/yperrot**(1./3.) 140 dhperrotdn=dhperrotduperrot*duperrotdyperrot*dyperrotdn 141 end if 142 143 etfw=hperrot*grho2_updn(ipts,1)*rho_inv*rho_inv 144 vtfw=-etfw + rho/hperrot*dhperrotdn*etfw 145 146 if(yperrot<=yperrot0)then 147 exci(ipts) = exci(ipts) + etfw + 1.5_dp*yperrot*dhperrotdyperrot*grho2_updn(ipts,1)*rho_inv*rho_inv 148 else 149 exci(ipts) = exci(ipts) + etfw + uperrot*dhperrotduperrot*grho2_updn(ipts,1)*rho_inv*rho_inv 150 end if 151 vxci(ipts,1) = vxci(ipts,1) + vtfw 152 if (has_fxc) fxci(ipts) = fxci(ipts) + etfw 153 if (has_dvxcdgr) dvxcdgr(ipts,1)= dvxcdgr(ipts,1)+two*hperrot*rho_inv 154 155 end do 156 157 ABI_DEALLOCATE(rho_updnm1_3) 158 159 !DBG_EXIT('COLL') 160 161 end subroutine xctfw