The Elk Code
xc_pzca.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 !BOP
7 ! !ROUTINE: xc_pzca
8 ! !INTERFACE:
9 subroutine xc_pzca(n,rho,ex,ec,vx,vc)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! n : number of density points (in,integer)
12 ! rho : charge density (in,real(n))
13 ! ex : exchange energy density (out,real(n))
14 ! ec : correlation energy density (out,real(n))
15 ! vx : exchange potential (out,real(n))
16 ! vc : correlation potential (out,real(n))
17 ! !DESCRIPTION:
18 ! Spin-unpolarised exchange-correlation potential and energy of the
19 ! Perdew-Zunger parameterisation of Ceperley-Alder electron gas: {\it Phys.
20 ! Rev. B} {\bf 23}, 5048 (1981) and {\it Phys. Rev. Lett.} {\bf 45}, 566
21 ! (1980).
22 !
23 ! !REVISION HISTORY:
24 ! Created October 2002 (JKD)
25 !EOP
26 !BOC
27 implicit none
28 ! arguments
29 integer, intent(in) :: n
30 real(8), intent(in) :: rho(n)
31 real(8), intent(out) :: ex(n),ec(n),vx(n),vc(n)
32 ! local variables
33 integer i
34 real(8), parameter :: pi=3.1415926535897932385d0
35 real(8), parameter :: thrd=1.d0/3.d0, thrd2=2.d0/3.d0, thrd4=4.d0/3.d0
36 real(8), parameter :: g=-0.1423d0,b1=1.0529d0,b2=0.3334d0
37 real(8), parameter :: a=0.0311d0,b=-0.048d0,c=0.0020d0,d=-0.0116d0
38 real(8) p1,p2,r,rs,t1
39 if (n < 1) then
40  write(*,*)
41  write(*,'("Error(xc_pzca): n < 1 : ",I8)') n
42  write(*,*)
43  stop
44 end if
45 ! prefactors
46 t1=3.d0/(4.d0*pi)
47 p1=t1**thrd
48 p2=t1*(9.d0*pi/4.d0)**thrd
49 do i=1,n
50  r=rho(i)
51  if (r < 1.d-12) then
52  ex(i)=0.d0
53  ec(i)=0.d0
54  vx(i)=0.d0
55  vc(i)=0.d0
56  cycle
57  end if
58  rs=p1/r**thrd
59 ! exchange energy and potential
60  ex(i)=-p2/rs
61  vx(i)=thrd4*ex(i)
62 ! correlation energy and potential
63  if (rs >= 1.d0) then
64  t1=sqrt(rs)
65  ec(i)=g/(1.d0+b1*t1+b2*rs)
66  vc(i)=ec(i)*(1.d0+(7.d0/6.d0)*b1*t1+thrd4*b2*rs)/(1.d0+b1*t1+b2*rs)
67  else
68  t1=dlog(rs)
69  ec(i)=a*t1+b+c*rs*t1+d*rs
70  vc(i)=a*t1+(b-thrd*a)+thrd2*c*rs*t1+thrd*(2.d0*d-c)*rs
71  end if
72 end do
73 end subroutine
74 !EOC
75 
subroutine xc_pzca(n, rho, ex, ec, vx, vc)
Definition: xc_pzca.f90:10