The Elk Code
x_wc06.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2006 Zhigang Wu and R. E. Cohen.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 elemental subroutine x_wc06(rho,s,u,v,ex,vx)
7 implicit none
8 ! arguments
9 real(8), intent(in) :: rho,s,u,v
10 real(8), intent(out) :: ex,vx
11 ! local variables
12 real(8), parameter :: ax=-0.7385587663820224059d0
13 real(8), parameter :: mu=0.2195149727645171d0
14 real(8), parameter :: kappa=0.804d0
15 real(8), parameter :: b=10.d0/81.d0
16 real(8), parameter :: c=0.00793746933516d0
17 real(8), parameter :: thrd=1.d0/3.d0
18 real(8), parameter :: thrd4=4.d0/3.d0
19 real(8) dmu,exu
20 real(8) s2,s4,es2,x,p0,fxwc
21 real(8) fs,fss,t0,t1,t2,t3
22 ! lda exchange energy density
23 exu=ax*rho**thrd
24 s2=s**2
25 s4=s2**2
26 es2=exp(-s2)
27 t0=1.d0+c*s4
28 dmu=mu-b
29 x=b*s2+dmu*s2*es2+log(t0)
30 p0=1.d0+x/kappa
31 ! WC enhancement factor
32 fxwc=1.d0+kappa-kappa/p0
33 ! exchange energy density
34 ex=exu*fxwc
35 t1=b+dmu*(1.d0-s2)*es2+2.d0*c*s2/t0
36 t2=dmu*s*(s2-2.d0)*es2+2.d0*c/t0-4.d0*(c**2)*s4/(t0**2)
37 t3=1.d0/(p0**2)
38 fs=2.d0*t1*t3
39 fss=t3*(4.d0*t2-8.d0*s*(t1**2)/(kappa*p0))
40 ! exchange potential
41 vx=exu*(thrd4*fxwc-(u-thrd4*s2*s)*fss-v*fs)
42 end subroutine
43 
elemental subroutine x_wc06(rho, s, u, v, ex, vx)
Definition: x_wc06.f90:7