TABLE OF CONTENTS
ABINIT/m_xcpositron [ Modules ]
NAME
m_xcpositron
FUNCTION
Compute electron-positron correlation potentials and energy density.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (GJ,MT) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
PARENTS
CHILDREN
SOURCE
21 #if defined HAVE_CONFIG_H 22 #include "config.h" 23 #endif 24 25 #include "abi_common.h" 26 27 module m_xcpositron 28 29 use defs_basis 30 use m_errors 31 use m_abicore 32 33 use m_numeric_tools, only : invcb 34 35 implicit none 36 37 private
ABINIT/xcpositron [ Functions ]
NAME
xcpositron
FUNCTION
Compute electron-positron correlation potentials and energy density. Used electron-positron correlation functional is controlled by ixcpositron argument. Returns Fxc, Vxc_pos, Vxc_el from input rhor_pos and rhor_el for positron and electrons.
INPUTS
grhoe2(ngr)=square of the gradient of electronic density rhoe (needed for GGA) ixcpositron=type of electron-positron correlation functional: 1: LDA zero positron density limit parametrized by Arponen & Pajanne and provided by Boronski & Nieminen [1,2] 11: LDA zero positron density limit parametrized by Arponen & Pajanne and fitted by Sterne & Kaiser [1,3] 2: LDA electron-positron correlation provided by Puska, Seitsonen, and Nieminen [1,4] 3: GGA zero positron density limit parametrized by Arponen & Pajanne and provided by Boronski & Nieminen [1,2,5] 31: GGA zero positron density limit parametrized by Arponen & Pajanne and fitted by Sterne & Kaiser [1,3,5] See references below ngr=size of grho2 array (0 if LDA, npt if GGA) npt=number of real space points on which density is provided posdensity0_limit=True if we are in the zero positron density limit rhoer(npt)=electron density (bohr^-3) rhopr(npt)=positron density (bohr^-3)
OUTPUT
fnxc(npt)=correlation energy per unit volume fxc vxce(npt)=correlation potential for electron dfxc/drhoe (hartree) vxcp(npt)=correlation potential for positron dfxc/drhop (hartree) vxcegr(ngr)= 1/|gradRhoe| dfxc/d|gradRhoe| (empty if LDA, i.e. ngr=0) Optional outputs: dvxce(npt)=partial second derivatives of the xc energy wr to the electronic density dvxce(:)=dVxce/dRhoe dvxcp(npt)=partial second derivatives of the xc energy wr to the positronic density dvxcp(:)=dVxcp/drhop
NOTES
References for electron-positron correlation functionals: [1] J. Arponen and E. Pajanne, Ann. Phys. (N.Y.) 121, 343 (1979) [[cite:Arponen1979a]]. [2] E. Boronski and R.M. Nieminen, Phys. Rev. B 34, 3820 (1986) [[cite:Boronski1986]]. [3] P.A. Sterne and J.H. Kaiser, Phys. Rev. B 43, 13892 (1991) [[cite:Sterne1991]]. [4] M.J. Puska, A.P. Seitsonen and R.M. Nieminen, Phys. Rev. B 52, 10947 (1994) [[cite:Puska1994]]. [5] B. Barbiellini, M.J. Puska, T. Torsti and R.M.Nieminen, Phys. Rev. B 51, 7341 (1995) [[cite:Barbiellini1995]]
PARENTS
m_pawxc,rhohxcpositron,rhotoxc
CHILDREN
invcb
SOURCE
102 subroutine xcpositron(fnxc,grhoe2,ixcpositron,ngr,npt,posdensity0_limit,rhoer,rhopr,vxce,vxcegr,vxcp,& 103 & dvxce,dvxcp) ! optional arguments 104 105 106 !This section has been created automatically by the script Abilint (TD). 107 !Do not modify the following lines by hand. 108 #undef ABI_FUNC 109 #define ABI_FUNC 'xcpositron' 110 !End of the abilint section 111 112 implicit none 113 114 !Arguments ------------------------------------ 115 !scalars 116 integer,intent(in) :: ixcpositron,ngr,npt 117 logical,intent(in) :: posdensity0_limit 118 !arrays 119 real(dp),intent(in) :: grhoe2(ngr),rhoer(npt),rhopr(npt) 120 real(dp),intent(out) :: fnxc(npt),vxce(npt),vxcegr(ngr),vxcp(npt) 121 real(dp),intent(out),optional :: dvxce(npt),dvxcp(npt) 122 123 !Local variables------------------------------- 124 !scalars 125 integer,parameter :: idebug=0 126 integer :: ipt 127 logical :: gga,need_dvxce,need_dvxcp 128 real(dp),parameter :: alpha_gga=0.22_dp 129 real(dp),parameter :: ap_a1=-1.56_dp,ap_b1=0.051_dp,ap_c1=-0.081_dp,ap_d1=1.14_dp 130 real(dp),parameter :: ap_a2=-0.92305_dp,ap_b2=-0.05459_dp 131 real(dp),parameter :: ap_a3=-0.6298_dp,ap_b3=-13.15111_dp,ap_c3=2.8655_dp 132 real(dp),parameter :: ap_a4=-179856.2768_dp,ap_b4=186.4207_dp,ap_c4=-0.524_dp 133 real(dp),parameter :: ap_psn_limit=0.7_dp 134 real(dp),parameter :: ap_psn_1=0.9_dp*ap_psn_limit,ap_psn_2=1.1_dp*ap_psn_limit 135 real(dp),parameter :: fpi3=third*four_pi 136 real(dp),parameter :: psn_aa=69.7029_dp,psn_ba=-107.4927_dp,psn_bb=141.8458_dp 137 real(dp),parameter :: psn_ca=23.7182_dp,psn_cb=-33.6472_dp ,psn_cc=5.21152_dp 138 real(dp),parameter :: sk_a=-1.56_dp,sk_b=0.1324_dp,sk_c=-4.092_dp,sk_d=51.96_dp,sk_e=0.7207_dp 139 real(dp),parameter :: rsfac=0.6203504908994000_dp 140 real(dp) :: arse,brse,crse,darse,dbrse,dcrse,d2arse,d2brse,d2crse 141 real(dp) :: d2eps,deps,dexc,dexcdg,dexc_p,d2expgga,dexpgga,d2invrs,dinvrs,d2kf,dkf,d2nqtf2,dnqtf2 142 real(dp) :: drse,drsp,d2exc,d2exc_p,d2rse,d2rsp,d2sqr,dsqr 143 real(dp) :: eexp,eps,exc,exc_p,expgga,invf,dinvf,d2invf,invrhoe,invrhop,invrs,invsqr 144 real(dp) :: kf,logrs,nqtf2,opr2,ratio_ap,ratio_psn,rhoe,rhop,rse,rsp,sqr 145 character(len=500) :: msg 146 !arrays 147 real(dp),allocatable :: rsepts(:),rsppts(:) 148 149 ! ************************************************************************* 150 151 gga=(ngr==npt) 152 need_dvxce=present(dvxce) 153 need_dvxcp=present(dvxcp) 154 155 if (gga.and.ixcpositron==2) then 156 msg = 'xcpositron: GGA not yet implemented for ixcpositron=2 !' 157 MSG_ERROR(msg) 158 end if 159 if (posdensity0_limit.and.ixcpositron==2) then 160 msg = 'xcpositron: ixcpositron=2 cannot be treated in the zero positron density limit !' 161 MSG_ERROR(msg) 162 end if 163 if (abs(ixcpositron)/=1.and.ixcpositron/=11.and.ixcpositron/=2.and.ixcpositron/=3.and.ixcpositron/=31) then 164 msg = 'xcpositron: unknown electron-positron correlation !' 165 MSG_ERROR(msg) 166 end if 167 168 !Compute density radii for rhor_el, rhor_pos 169 ABI_ALLOCATE(rsepts,(npt)) 170 call invcb(rhoer(:),rsepts,npt) 171 rsepts(:)=rsfac*rsepts(:) 172 if (ixcpositron==2) then 173 ABI_ALLOCATE(rsppts,(npt)) 174 call invcb(rhopr(:),rsppts,npt) 175 rsppts(:)=rsfac*rsppts(:) 176 end if 177 178 !Loop over grid points 179 !---------------------------------------------------- 180 do ipt=1,npt 181 182 rhoe=rhoer(ipt) 183 rhop=rhopr(ipt) 184 exc=zero;dexc=zero;d2exc=zero;dexcdg=zero 185 186 rse=rsepts(ipt) 187 invrhoe=one/rhoe 188 drse=-third*rse*invrhoe 189 if (need_dvxce) d2rse= four/nine*rse*invrhoe**2 190 191 ! Arponen & Pajane parametrization for electron 192 if (ixcpositron/=11.and.ixcpositron/=31) then 193 if (rse<0.302_dp) then 194 invrs=one/rse;invsqr=sqrt(invrs);logrs=log(rse) 195 exc =ap_a1*invsqr+(ap_b1*logrs+ap_c1)*logrs+ap_d1 196 dexc=drse*invrs*(-half*ap_a1*invsqr+two*ap_b1*logrs+ap_c1) 197 if (need_dvxce) d2exc=(d2rse/drse-drse*invrs)*dexc+drse**2*invrs**2*(quarter*ap_a1*invsqr+two*ap_b1) 198 else if (rse>=0.302_dp.and.rse<=0.56_dp) then 199 invrs=one/rse 200 exc =ap_a2+ap_b2*invrs**2 201 dexc=-drse*ap_b2*two*invrs**3 202 if (need_dvxce) d2exc=d2rse/drse*dexc+six*drse**2*ap_b2*invrs**4 203 else if (rse>0.56_dp.and.rse<=8.0_dp) then 204 invrs=one/(rse+2.5_dp) 205 dinvrs=-drse*invrs**2 206 ! jmb : d2rse initialized only if need_dvxce = .True. 207 if (need_dvxce) d2invrs=-d2rse*invrs**2-two*invrs*drse**2 208 exc =ap_a3+ap_b3*invrs**2+ap_c3*invrs 209 dexc=two*ap_b3*invrs*dinvrs+ap_c3*dinvrs 210 if (need_dvxce) d2exc=two*ap_b3*dinvrs**2+(two*ap_b3*invrs+ap_c3)*d2invrs 211 else 212 exc =ap_a4*rhoe**2+ap_b4*rhoe+ap_c4 213 dexc =two*ap_a4*rhoe+ap_b4 214 if (need_dvxce) d2exc=two*ap_a4 215 end if 216 217 ! Sterne & Kaiser parametrization for electron 218 else 219 eexp=exp(-(rse+sk_c)**2/sk_d) 220 opr2=(one+rse**2) 221 arse=atan(rse) 222 exc = sk_a/sqrt(arse)+sk_b*eexp+sk_e 223 dexc= -(two*sk_b*eexp*(sk_c+rse)/sk_d + sk_a/(two*opr2*sqrt(arse)**3))*drse 224 if (need_dvxce) d2exc=-(two*sk_b*eexp*(sk_c+rse)/sk_d + sk_a/(two*opr2*arse**1.5_dp))*d2rse & 225 & +(two*sk_b*eexp*(two*sk_c**2-sk_d+four*sk_c*rse+two*rse**2)/sk_d**2 & 226 & +sk_a*(three+four*rse*arse)/(four*opr2**2*sqrt(arse)**5))*drse**2 227 end if 228 229 ! Puska, Seitsonen and Nieminen parametrization for positron 230 if (ixcpositron==2.and.rse>=ap_psn_1) then 231 rsp=rsppts(ipt) 232 invrhop=one/rhop 233 drsp=-third*rsp*invrhop 234 if (need_dvxcp) d2rsp= four/nine*rsp*invrhop**2 235 exc_p=zero;dexc_p=zero;d2exc_p=zero 236 if (rsp<0.302_dp) then 237 invrs=one/rsp;invsqr=sqrt(invrs);logrs=log(rsp) 238 exc_p =ap_a1*invsqr+(ap_b1*logrs+ap_c1)*logrs+ap_d1 239 dexc_p=drsp*invrs*(-half*ap_a1*invsqr+two*ap_b1*logrs+ap_c1) 240 if (need_dvxcp) d2exc_p=(d2rsp/drsp-drsp*invrs)*dexc_p+drsp**2*invrs**2*(quarter*ap_a1*invsqr+two*ap_b1) 241 else if (rsp>=0.302_dp.and.rsp<=0.56_dp) then 242 invrs=one/rsp 243 exc_p =ap_a2+ap_b2*invrs**2 244 dexc_p=-drsp*ap_b2*two*invrs**3 245 if (need_dvxcp) d2exc_p=d2rsp/drsp*dexc_p+six*drsp**2*ap_b2*invrs**4 246 else if (rsp>0.56_dp.and.rsp<=8.0_dp) then 247 invrs=one/(rsp+2.5_dp) 248 dinvrs=-drsp*invrs**2 249 ! jmb : d2rsp initialized only if need_dvxcp = .True.* 250 if (need_dvxcp) d2invrs=-d2rsp*invrs**2-two*invrs*drsp**2 251 exc_p =ap_a3+ap_b3*invrs**2+ap_c3*invrs 252 dexc_p=two*ap_b3*invrs*dinvrs+ap_c3*dinvrs 253 if (need_dvxcp) d2exc_p=two*ap_b3*dinvrs**2+(two*ap_b3*invrs+ap_c3)*d2invrs 254 else 255 exc_p =ap_a4*rhop**2+ap_b4*rhop+ap_c4 256 dexc_p =two*ap_a4*rhop+ap_b4 257 if (need_dvxcp) d2exc_p=two*ap_a4 258 end if 259 end if 260 261 ! GGA correction 262 if (gga) then 263 kf=(three*pi*pi*rhoe)**third 264 nqtf2=(rhoe*sqrt(four*kf/pi))**2 265 eps=grhoe2(ipt)/nqtf2 266 if (eps<zero) then 267 MSG_ERROR('xcpositron: problem, negative GGA espilon !') 268 end if 269 expgga=exp(-alpha_gga*eps*third) 270 271 dkf=pi*pi/(sqrt(three*pi*pi*rhoe)**third) 272 d2kf=-two*pi*pi*pi*pi*(three*pi*pi*rhoe)**(-5.0_dp/3.0_dp) 273 sqr=sqrt(four*kf/pi) 274 dsqr=(four*dkf/pi)/(two*sqr) 275 d2sqr=two/(pi*sqr*dkf)*(d2kf*sqr-dsqr*dkf) 276 nqtf2=(rhoe*sqr)**two 277 dnqtf2=two*(sqr+rhoe*dsqr)*rhoe*sqr 278 d2nqtf2=two*(rhoe*sqr*(two*dsqr+rhoe*d2sqr) & 279 & +sqr*(sqr+rhoe*dsqr) & 280 & +rhoe*(sqr+rhoe*dsqr) ) 281 deps=-grhoe2(ipt)*dnqtf2/(nqtf2**two) 282 d2eps=-grhoe2(ipt)/(nqtf2*nqtf2*dnqtf2)*(d2nqtf2*nqtf2*nqtf2-two*nqtf2*dnqtf2*dnqtf2) 283 dexpgga=-alpha_gga*third*deps*expgga 284 d2expgga=-alpha_gga*third*(d2eps*expgga+deps*dexpgga) 285 286 exc = exc *expgga 287 dexc=(dexc*expgga+exc*dexpgga) 288 if (need_dvxce) d2exc=d2exc*expgga+two*dexc*dexpgga+exc*d2expgga 289 if (abs(grhoe2(ipt))<1.e24_dp) dexcdg=-exc*alpha_gga*two_thirds/nqtf2 290 end if 291 292 ! Computation of XC energy, potentials and kernels 293 ! Zero positron density limit 294 if (ixcpositron/=2.or.rse<ap_psn_1) then 295 fnxc(ipt)=rhop*exc 296 vxce(ipt)=rhop*dexc 297 vxcp(ipt)=exc 298 if (need_dvxce) dvxce(ipt)=rhop*d2exc 299 if (need_dvxcp) dvxcp(ipt)=zero 300 if (gga) vxcegr(ipt)=rhop*dexcdg 301 else 302 ! Puska, Seitsonen and Nieminen functional 303 arse=psn_aa+psn_ba*rse+psn_ca*rse**2 304 brse=psn_ba+psn_bb*rse+psn_cb*rse**2 305 crse=psn_ca+psn_cb*rse+psn_cc*rse**2 306 darse=(psn_ba+two*psn_ca*rse)*drse 307 dbrse=(psn_bb+two*psn_cb*rse)*drse 308 dcrse=(psn_cb+two*psn_cc*rse)*drse 309 invf=arse+brse*rsp+crse*rsp**2+invrhop/exc+invrhoe/exc_p 310 fnxc(ipt)=one/invf 311 dinvf=darse+dbrse*rsp+dcrse*rsp**2-invrhop*dexc/exc**2-invrhoe**2/exc_p 312 vxce(ipt)=-dinvf/invf**2 313 if (need_dvxce) then 314 d2arse=darse*d2rse/drse+two*psn_ca*drse**2 315 d2brse=dbrse*d2rse/drse+two*psn_cb*drse**2 316 d2crse=dcrse*d2rse/drse+two*psn_cc*drse**2 317 d2invf=d2arse+d2brse*rsp+d2crse*rsp**2 & 318 & +invrhop*(two*dexc**2/exc-d2exc)/exc**2+two*invrhoe**3/exc_p 319 dvxce(ipt)=(two*dinvf**2/invf-d2invf)/invf**2 320 end if 321 dinvf=(brse+two*crse*rsp)*drsp-invrhop**2/exc-invrhoe*dexc_p/exc_p**2 322 vxcp(ipt)=-dinvf/invf**2 323 if (need_dvxcp) then 324 d2invf=two*crse*drsp+(brse+two*crse*rsp)*d2rsp & 325 & +two*invrhop**3/exc+invrhoe*(two*dexc_p**2/exc_p-d2exc_p)/exc_p**2 326 dvxcp(ipt)=(two*dinvf**2/invf-d2invf)/invf**2 327 end if 328 ! For small rse, use pure Arponen/Pajanne functional 329 ! Around the limit (rse=0.7, see PSN paper), switch smoothly from PSN to AP 330 if (rse>=ap_psn_1.and.rse<=ap_psn_2) then 331 ratio_psn=(rse-ap_psn_1)/(ap_psn_2-ap_psn_1);ratio_ap=one-ratio_psn 332 fnxc(ipt)=ratio_psn*fnxc(ipt)+ratio_ap*rhop*exc 333 vxce(ipt)=ratio_psn*vxce(ipt)+ratio_ap*rhop*dexc 334 vxcp(ipt)=ratio_psn*vxcp(ipt)+ratio_ap*exc 335 if (need_dvxce) dvxce(ipt)=ratio_psn*dvxce(ipt)+ratio_ap*rhop*d2exc 336 if (need_dvxcp) dvxcp(ipt)=ratio_psn*dvxcp(ipt) 337 end if 338 end if 339 340 ! Debug statements: use polynomial functionals 341 if (idebug>0) then 342 if (idebug==4) then ! order 4 343 fnxc(ipt)=tol3*((rhop**4+rhoe**4)/12._dp+(rhop**3*rhoe+rhop*rhoe**3)/3._dp+rhop**2*rhoe**2) 344 vxce(ipt)=tol3*((rhop**3*rhoe+rhop*rhoe**3)/3._dp+rhop**2*rhoe+rhop*rhoe**2) 345 vxcp(ipt)=tol3*((rhop**3*rhoe+rhop*rhoe**3)/3._dp+rhop**2*rhoe+rhop*rhoe**2) 346 if (need_dvxce) dvxce(ipt)=tol3*(rhop**3/3._dp+rhop*rhoe**2+rhop**2+two*rhop*rhoe) 347 if (need_dvxcp) dvxcp(ipt)=tol3*(rhoe**3/3._dp+rhoe*rhop**2+rhoe**2+two*rhop*rhoe) 348 end if 349 if (idebug==3) then ! order 3 350 fnxc(ipt)=tol3*((rhop**3+rhoe**3)*third+rhop**2*rhoe+rhop*rhoe**2) 351 vxce(ipt)=tol3*(rhop+rhoe)**2 352 vxcp(ipt)=tol3*(rhop+rhoe)**2 353 if (need_dvxce) dvxce(ipt)=tol3*two*rhoe 354 if (need_dvxcp) dvxcp(ipt)=tol3*two*rhop 355 end if 356 if (idebug==2) then ! order 2 357 fnxc(ipt)=tol3*(rhop+rhoe)**2 358 vxce(ipt)=tol3*two*(rhop+rhoe) 359 vxcp(ipt)=tol3*two*(rhop+rhoe) 360 if (need_dvxce) dvxce(ipt)=tol3*two 361 if (need_dvxcp) dvxcp(ipt)=tol3*two 362 end if 363 if (idebug==1) then ! order 1 364 fnxc(ipt)=tol3*(rhop+rhoe) 365 vxce(ipt)=tol3 366 vxcp(ipt)=tol3 367 if (need_dvxce) dvxce(ipt)=zero 368 if (need_dvxcp) dvxcp(ipt)=zero 369 end if 370 end if 371 372 end do ! ipt 373 374 ABI_DEALLOCATE(rsepts) 375 if (ixcpositron==2) then 376 ABI_DEALLOCATE(rsppts) 377 end if 378 379 !Convert everything in Hartree units 380 fnxc(:)=half*fnxc(:) 381 vxce(:)=half*vxce(:) 382 vxcp(:)=half*vxcp(:) 383 if (need_dvxce) dvxce(:)=half*dvxce(:) 384 if (need_dvxcp) dvxcp(:)=half*dvxcp(:) 385 if (gga) vxcegr(:)=half*vxcegr(:) 386 387 end subroutine xcpositron