TABLE OF CONTENTS
ABINIT/xclb [ Functions ]
NAME
xclb
FUNCTION
Computes the GGA like part (vx_lb) of the Leeuwen-Baerends exchange-correlation potential (vxc_lb) and adds it to the lda exchange-correlation potential (vxc_lda) which must be provided as input, vxci <-- vxc_lb =: vxc_lda + vx_lb [R van Leeuwen and EJ Baerends, Phys Rev A 49, 2421 (1994)] With respect to spin, the van Leeuwen-Baerends potential is "exchange-like" : separate contributions from spin up and spin down.
COPYRIGHT
Copyright (C) 2000-2018 ABINIT group (MF,LG) 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.txt .
INPUTS
npts= number of points to be computed nspden=1 for unpolarized, 2 for spin-polarized grho2_updn(npts,2*nspden-1)=square of the gradient of the spin-up, and, if nspden==2, spin-down, and total density (Hartree/Bohr**2) rho_updn(npts,nspden)=spin-up and spin-down density (Hartree/bohr**3)
OUTPUT
(see side effects)
SIDE EFFECTS
Input/Output: vxci(npts,nspden)=input xc potential to which Leeuwen-Baerends correction is added at output.
PARENTS
drivexc
CHILDREN
SOURCE
48 #if defined HAVE_CONFIG_H 49 #include "config.h" 50 #endif 51 52 #include "abi_common.h" 53 54 55 subroutine xclb(grho2_updn,npts,nspden,rho_updn,vxci) 56 57 use defs_basis 58 use m_profiling_abi 59 60 !This section has been created automatically by the script Abilint (TD). 61 !Do not modify the following lines by hand. 62 #undef ABI_FUNC 63 #define ABI_FUNC 'xclb' 64 !End of the abilint section 65 66 implicit none 67 68 !Arguments ------------------------------------ 69 !scalars 70 integer,intent(in) :: npts,nspden 71 !arrays 72 real(dp),intent(in) :: grho2_updn(npts,2*nspden-1),rho_updn(npts,nspden) 73 real(dp),intent(inout) :: vxci(npts,nspden) 74 75 !Local variables------------------------------- 76 !scalars 77 integer :: ipts,ispden 78 real(dp),parameter :: beta=0.05_dp 79 real(dp) :: density,density_gradient,density_t13,s_g_sq,scaled_gradient 80 real(dp) :: scaling_factor,vx_lb 81 82 ! ************************************************************************* 83 84 !DEBUG 85 !write(std_out,*) ' %xclb: enter' 86 !ENDDEBUG 87 88 !scale the spin densities for evaluating spin up or down exchange 89 scaling_factor=one 90 if(nspden == 2) scaling_factor=two 91 92 do ispden=1,nspden 93 94 do ipts=1,npts 95 96 density= scaling_factor * rho_updn(ipts,ispden) 97 density_gradient= scaling_factor * sqrt(grho2_updn(ipts,ispden)) 98 99 density_t13= density**third 100 scaled_gradient= density_gradient/max(density*density_t13,1.e-12_dp) 101 102 s_g_sq= scaled_gradient*scaled_gradient 103 104 vx_lb= -beta*density_t13 * s_g_sq/ & 105 & (one+3.d0*beta* scaled_gradient*log(scaled_gradient+sqrt(one+s_g_sq*s_g_sq))) 106 107 vxci(ipts,ispden)=vxci(ipts,ispden)+vx_lb 108 end do 109 110 end do 111 112 end subroutine xclb