The Elk Code
xc_am05.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2004, 2005 Rickard Armiento
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_am05
8 ! !INTERFACE:
9 subroutine xc_am05(n,rho,grho,g2rho,g3rho,ex,ec,vx,vc)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! n : number of density points (in,integer)
12 ! rho : charge density (in,real(n))
13 ! grho : |grad rho| (in,real(n))
14 ! g2rho : grad^2 rho (in,real(n))
15 ! g3rho : (grad rho).(grad |grad rho|) (in,real(n))
16 ! ex : exchange energy density (out,real(n))
17 ! ec : correlation energy density (out,real(n))
18 ! vx : spin-unpolarised exchange potential (out,real(n))
19 ! vc : spin-unpolarised correlation potential (out,real(n))
20 ! !DESCRIPTION:
21 ! Spin-unpolarised exchange-correlation potential and energy functional of
22 ! R. Armiento and A. E. Mattsson, {\it Phys. Rev. B} {\bf 72}, 085108 (2005).
23 !
24 ! !REVISION HISTORY:
25 ! Created April 2005 (RAR); based on xc_pbe
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 integer, intent(in) :: n
31 real(8), intent(in) :: rho(n),grho(n),g2rho(n),g3rho(n)
32 real(8), intent(out) :: ex(n),ec(n),vx(n),vc(n)
33 ! local variables
34 integer i
35 real(8), parameter :: pi=3.1415926535897932385d0
36 real(8) r,kf,s,v,u
37 real(8) grho_,g2rho_,g3rho_
38 do i=1,n
39  r=rho(i)
40  if (r > 1.d-12) then
41  grho_=grho(i)
42  g2rho_=g2rho(i)
43  g3rho_=g3rho(i)
44 ! exchange energy density and potential
45  kf=(r*3.d0*pi**2)**(1.d0/3.d0)
46  s=grho_/(2.d0*kf*r)
47  v=g2rho_/(r*(2.d0*kf)**2)
48  u=g3rho_/((r**2)*(2.d0*kf)**3)
49  call xc_am05_point(r,s,u,v,ex(i),ec(i),vx(i),vc(i),1)
50  else
51  ex(i)=0.d0
52  ec(i)=0.d0
53  vx(i)=0.d0
54  vc(i)=0.d0
55  end if
56 end do
57 end subroutine
58 !EOC
59 
60 !BOP
61 ! !ROUTINE: xc_am05_point
62 ! !INTERFACE:
63 subroutine xc_am05_point(rho,s,u,v,ex,ec,vx,vc,pot)
64 ! !INPUT/OUTPUT PARAMETERS:
65 ! rho : electron density (in,real)
66 ! s : gradient of n / (2 kF n)
67 ! u : grad n * grad | grad n | / (n**2 (2 kF)**3)
68 ! v : laplacian of density / (n**2 (2.d0*kf)**3)
69 ! ex : exchange energy density (out,real)
70 ! ec : correlation energy density (out,real)
71 ! vx : spin-unpolarised exchange potential (out,real)
72 ! vc : spin-unpolarised correlation potential (out,real)
73 ! !DESCRIPTION:
74 ! Calculate the spin-unpolarised exchange-correlation potential and energy for
75 ! the Armiento-Mattsson 05 functional for a single point.
76 !
77 ! !REVISION HISTORY:
78 ! Created April 2005 (RAR)
79 !EOP
80 !BOC
81 implicit none
82 ! arguments
83 real(8), intent(in) :: rho, s, u, v
84 integer, intent(in) :: pot
85 real(8), intent(out) :: ex, ec, vx, vc
86 ! constants
87 real(8), parameter :: pi=3.1415926535897932385d0
88 real(8), parameter :: g = 0.8098d0
89 real(8), parameter :: a = 2.804d0
90 real(8), parameter :: c = 0.7168d0
91 ! local variables
92 real(8) s2,exlda, vxlda, eclda, vclda, X, Xs, Xss
93 real(8) F, Fs, Fss, Hx, Hxs, Hxss, Hc, Hcs, Hcss
94 real(8) zb, zbs, zbss, w
95 real(8) n0b, n0bs, n0bss
96 real(8) ln0b, ln0bs, ln0bss
97 real(8) zbb, zbbc, zbbs, zbbss
98 real(8) fxb, fxbs, fxbss
99 ! cutoff
100 if((rho <= 1.d-16)) then
101  ex = 0.0d0
102  ec = 0.0d0
103  vx = 0.0d0
104  vc = 0.0d0
105  return
106 endif
107 s2 = s**2
108 ! LDA correlation
109 call xc_am05_ldapwc(rho,eclda,vclda)
110 ! LDA exchange
111 call xc_am05_ldax(rho,exlda,vxlda)
112 !------------------!
113 ! exchange !
114 !------------------!
115 ! interpolation index
116 x = 1.0d0 - a*s2/(1.0d0 + a*s2)
117 ! Airy LAA refinement function
118 call xc_am05_labertw(s**(3.0d0/2.0d0)/sqrt(24.0d0),w)
119 zb = (3.0d0/2.0d0*w)**(2.0d0/3.0d0)
120 n0b = w/(2.0d0*pi**2*s**3)
121 ln0b = -3.0d0/(2.0d0*pi)*(3.0d0*pi**2*n0b)**(1.0d0/3.0d0)
122 zbbc = ((4.0d0/3.0d0)**(1.0d0/3.0d0)*2.0d0*pi/3.0d0)**4
123 zbb = (zbbc*zb**2 + zb**4)**(1.0d0/4.0d0)
124 fxb = -1.0d0/(ln0b*2.0d0*zbb)
125 f = (c*s2 + 1.0d0)/(c*s2/fxb + 1.0d0)
126 ! exchange refinement function
127 hx = x + (1.0d0 - x)*f
128 ! exchange energy per particle, Ex = Integrate[n*ex]
129 ex = exlda*hx
130 !---------------------!
131 ! correlation !
132 !---------------------!
133 ! correlation refinement function
134 hc = x + g*(1.0d0 - x)
135 ! correlation energy per particle, Ec = Integrate[rho*ec]
136 ec = eclda*hc
137 if (pot == 0) return
138 !----------------------------!
139 ! exchange potential !
140 !----------------------------!
141 ! interpolation index derivatives, dX/ds
142 xs = -2.0d0*a*s/(1.0d0 + a*s2)**2
143 xss = 2.0d0*a*(3.0d0*a*s2-1.0d0)/(1.0d0+a*s2)**3
144 ! airy LAA refinement function derivatives, dF/ds
145 zbs = zb/(s + s*w)
146 zbss = - zb*w*(5.0d0+2.0d0*w)/(2.0d0*s2*(1.0d0+w)**3)
147 n0bs = sqrt(zb)*(-2.0d0*zb+s*zbs)/(2.0d0*pi**2*s2**2)
148 n0bss = (16.0d0*zb**2+s**2*zbs**2+2.0d0*s*zb*(-6.0d0* &
149  zbs+s*zbss))/(4.0d0*pi**2*s**5*sqrt(zb))
150 ln0bs = -(3.0d0/pi)**(1.0d0/3.0d0)*n0bs/ &
151  (2.0d0*n0b**(2.0d0/3.0d0))
152 ln0bss = (2.0d0*n0bs**2-3.0d0*n0b*n0bss)/(2.0d0* &
153  3.0d0**(2.0d0/3.0d0)*pi**(1.0d0/3.0d0)*n0b**(5.0d0/3.0d0))
154 zbbs = zb*(zbbc+2*zb**2)*zbs/ &
155  (2.0d0*(zb**2*(zbbc+zb**2))**(3.0d0/4.0d0))
156 zbbss = zb**2*(-zbbc*(zbbc-2.0d0*zb**2)*zbs**2+ &
157  2.0d0*zb*(zbbc+zb**2)*(zbbc+2.0d0*zb**2)*zbss)/ &
158  (4.0d0*(zb**2*(zbbc+zb**2))**(7.0d0/4.0d0))
159 fxbs = (zbb*ln0bs+ln0b*zbbs)/(2.0d0*ln0b**2*zbb**2)
160 fxbss = (-2.0d0*ln0b**2*zbbs**2+zbb**2*(-2.0d0*ln0bs**2 + &
161  ln0b*ln0bss)+ln0b*zbb*(-2.0d0*ln0bs*zbbs+ln0b*zbbss))/ &
162  (2.0d0*ln0b**3*zbb**3)
163 fs = (c*s*(2.0d0*(fxb-1.0d0)*fxb + s*(1.0d0+c*s2)*fxbs))/ &
164  (c*s2 + fxb)**2
165 fss = (c*(-2.0d0*(3.0d0*c*s2-fxb)*(fxb-1.0d0)*fxb+ &
166  4.0d0*s*(-c*s2+fxb+2.0d0*c*s2*fxb)*fxbs - &
167  2.0d0*s2*(1.0d0+c*s2)*fxbs**2+s2*(1.0d0+c*s2)* &
168  (c*s2 + fxb)*fxbss))/(c*s2+fxb)**3
169 ! GGA refinement function derivatives, dF/ds
170 hxs = - (x - 1.0d0)*fs - (f - 1.0d0)*xs
171 hxss = - 2.0d0*fs*xs - (x - 1.0d0)*fss - (f - 1.0d0)*xss
172 ! vx formula for gradient dependent functional,
173 ! generalized form of Eq. (24) in PRB 33, 8800 (1986)
174 vx = vxlda*(hx - s*hxs) + &
175  exlda*((4.0d0/3.0d0*s-v/s)*hxs - &
176  (u-4.0d0/3.0d0*s**3)*(hxss/s-hxs/s2))
177 !-------------------------------!
178 ! correlation potential !
179 !-------------------------------!
180 ! correlation refinement function derivatives, dF/ds
181 hcs = xs - g*xs
182 hcss = xss - g*xss
183 ! vc formula for gradient dependent functional,
184 ! generalized form of Eq. (24) in Phys. Rev. B 33, 8800 (1986)
185 vc = vclda*(hc - s*hcs) + &
186  eclda*((4.0d0/3.0d0*s - v/s)*hcs - &
187  (u - 4.0d0/3.0d0*s**3)*(hcss/s - hcs/s2))
188 end subroutine
189 !EOC
190 
191 !BOP
192 ! !ROUTINE: xc_am05_ldax
193 ! !INTERFACE:
194 subroutine xc_am05_ldax(n,ex,vx)
195 ! !INPUT/OUTPUT PARAMETERS:
196 ! n : electron density (in,real)
197 ! ex : exchange energy per electron (out,real)
198 ! vx : exchange potential (out,real)
199 ! !DESCRIPTION:
200 ! Local density approximation exchange.
201 !
202 ! !REVISION HISTORY:
203 ! Created April 2005 (RAR)
204 !EOP
205 !BOC
206 implicit none
207 ! arguments
208 real(8), intent(in) :: n
209 real(8), intent(out) :: ex
210 real(8), intent(out) :: vx
211 ! constants
212 real(8), parameter :: pi=3.1415926535897932385d0
213 vx=-(3.d0*n/pi)**(1.d0/3.d0)
214 ex=(3.d0/4.d0)*vx
215 end subroutine
216 !EOC
217 
218 !BOP
219 ! !ROUTINE: xc_am05_ldapwc
220 ! !INTERFACE:
221 subroutine xc_am05_ldapwc(n,ec,vc)
222 ! !INPUT/OUTPUT PARAMETERS:
223 ! n : electron density (in,real)
224 ! ec : correlation energy per electron (out,real)
225 ! vc : correlation potential (out,real)
226 ! !DESCRIPTION:
227 ! Correlation energy and potential of the Perdew-Wang parameterisation of
228 ! the Ceperley-Alder electron gas {\it Phys. Rev. B} {\bf 45}, 13244 (1992)
229 ! and {\it Phys. Rev. Lett.} {\bf 45}, 566 (1980). This is a clean-room
230 ! implementation from paper.
231 !
232 ! !REVISION HISTORY:
233 ! Created April 2005 (RAR)
234 !EOP
235 !BOC
236 implicit none
237 ! arguments
238 real(8), intent(in) :: n
239 real(8), intent(out) :: ec
240 real(8), intent(out) :: vc
241 ! constants
242 real(8), parameter :: pi=3.1415926535897932385d0
243 real(8), parameter :: a01 = 0.21370d0
244 real(8), parameter :: b01 = 7.5957d0
245 real(8), parameter :: b02 = 3.5876d0
246 real(8), parameter :: b03 = 1.6382d0
247 real(8), parameter :: b04 = 0.49294d0
248 ! paper actually use this:
249 ! real(8), parameter (A0 = 0.031091d0)
250 ! but routines now "defacto standard" was distributed using:
251 real(8), parameter :: A0 = 0.0310907d0
252 ! local variables
253 real(8) rsq
254 real(8) Q0, Q1, Q1p, ecrs
255 rsq = (3.0d0/(4.0d0*pi*n))**(1.0d0/6.0d0)
256 ec = -2.0d0*a0*(1.0d0 + a01*rsq**2)* &
257  log(1.0d0 + 1.0d0/ &
258  (2.0d0*a0*rsq*(b01 + rsq*(b02 + rsq*(b03 + b04*rsq)))))
259 q0 = -2.0d0*a0*(1.0d0 + a01*rsq**2)
260 q1 = 2.0d0*a0*rsq*(b01 + rsq*(b02 + rsq*(b03 + b04*rsq)))
261 q1p = a0*(b01/rsq+2.0d0*b02+3.0d0*b03*rsq+4.0d0*b04*rsq**2)
262 ecrs = -2.0d0*a0*a01*log(1.0d0 + 1.0d0/q1)-q0*q1p/(q1**2+q1)
263 vc = ec - rsq**2/3.0d0*ecrs
264 end subroutine
265 !EOC
266 
267 !BOP
268 ! !ROUTINE: xc_am05_labertw
269 ! !INTERFACE:
270 subroutine xc_am05_labertw(z,val)
271 ! !INPUT/OUTPUT PARAMETERS:
272 ! z : function argument (in,real)
273 ! val : value of lambert W function of z (out,real)
274 ! !DESCRIPTION:
275 ! Lambert $W$-function using the method of Corless, Gonnet, Hare, Jeffrey and
276 ! Knuth, {\it Adv. Comp. Math.} {\bf 5}, 329 (1996). The approach is based
277 ! loosely on that in GNU Octave by N. N. Schraudolph, but this implementation
278 ! is for real values and the principal branch only.
279 !
280 ! !REVISION HISTORY:
281 ! Created April 2005 (RAR)
282 !EOP
283 !BOC
284 implicit none
285 ! arguments
286 real(8), intent(in) :: z
287 real(8), intent(out) :: val
288 ! local variables
289 real(8) e,t,p
290 integer i
291 ! if z too low, go with the first term of the power expansion, z
292 if (z < 1.d-20) then
293  val=z
294  return
295 end if
296 e=exp(1.d0)
297 ! inital guess
298 if (abs(z+1.d0/e) > 1.45d0) then
299 ! asymptotic expansion at 0 and Inf
300  val=log(z)
301  val=val-log(val)
302 else
303 ! series expansion about -1/e to first order
304  val=1.d0*sqrt(2.d0*e*z+2.d0)-1.d0
305 end if
306 ! find val through iteration
307 do i=1,10
308  p=exp(val)
309  t=val*p-z
310  if (val /= -1.d0) then
311  t=t/(p*(val+1.d0)-0.5d0*(val+2.d0)*t/(val+1.d0))
312  else
313  t=0.d0
314  end if
315  val=val-t
316  if (abs(t) < (2.48d0*1.d-14)*(1.d0+abs(val))) return
317 end do
318 ! this should never happen!
319 write(*,*)
320 write(*,'("Error(xc_am05_labertw): iteration limit reached")')
321 write(*,'(" Likely cause: improper numbers (INFs, NaNs) in density")')
322 write(*,*)
323 stop
324 end subroutine
325 !EOC
326 
subroutine xc_am05_ldax(n, ex, vx)
Definition: xc_am05.f90:195
subroutine xc_am05_labertw(z, val)
Definition: xc_am05.f90:271
subroutine xc_am05(n, rho, grho, g2rho, g3rho, ex, ec, vx, vc)
Definition: xc_am05.f90:10
subroutine xc_am05_ldapwc(n, ec, vc)
Definition: xc_am05.f90:222
subroutine xc_am05_point(rho, s, u, v, ex, ec, vx, vc, pot)
Definition: xc_am05.f90:64