The Elk Code
 
Loading...
Searching...
No Matches
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
6elemental subroutine x_wc06(rho,s,u,v,ex,vx)
7implicit none
8! arguments
9real(8), intent(in) :: rho,s,u,v
10real(8), intent(out) :: ex,vx
11! local variables
12real(8), parameter :: ax=-0.7385587663820224059d0
13real(8), parameter :: mu=0.2195149727645171d0
14real(8), parameter :: kappa=0.804d0
15real(8), parameter :: b=10.d0/81.d0
16real(8), parameter :: c=0.00793746933516d0
17real(8), parameter :: thrd=1.d0/3.d0
18real(8), parameter :: thrd4=4.d0/3.d0
19real(8) dmu,exu
20real(8) s2,s4,es2,x,p0,fxwc
21real(8) fs,fss,t0,t1,t2,t3
22! lda exchange energy density
23exu=ax*rho**thrd
24s2=s**2
25s4=s2**2
26es2=exp(-s2)
27t0=1.d0+c*s4
28dmu=mu-b
29x=b*s2+dmu*s2*es2+log(t0)
30p0=1.d0+x/kappa
31! WC enhancement factor
32fxwc=1.d0+kappa-kappa/p0
33! exchange energy density
34ex=exu*fxwc
35t1=b+dmu*(1.d0-s2)*es2+2.d0*c*s2/t0
36t2=dmu*s*(s2-2.d0)*es2+2.d0*c/t0-4.d0*(c**2)*s4/(t0**2)
37t3=1.d0/(p0**2)
38fs=2.d0*t1*t3
39fss=t3*(4.d0*t2-8.d0*s*(t1**2)/(kappa*p0))
40! exchange potential
41vx=exu*(thrd4*fxwc-(u-thrd4*s2*s)*fss-v*fs)
42end subroutine
43
elemental subroutine x_wc06(rho, s, u, v, ex, vx)
Definition x_wc06.f90:7