9subroutine xc_am05(n,rho,grho,g2rho,g3rho,ex,ec,vx,vc)
30integer,
intent(in) :: n
31real(8),
intent(in) :: rho(n),grho(n),g2rho(n),g3rho(n)
32real(8),
intent(out) :: ex(n),ec(n),vx(n),vc(n)
35real(8),
parameter :: pi=3.1415926535897932385d0
37real(8) grho_,g2rho_,g3rho_
45 kf=(r*3.d0*pi**2)**(1.d0/3.d0)
47 v=g2rho_/(r*(2.d0*kf)**2)
48 u=g3rho_/((r**2)*(2.d0*kf)**3)
83real(8),
intent(in) :: rho, s, u, v
84integer,
intent(in) :: pot
85real(8),
intent(out) :: ex, ec, vx, vc
87real(8),
parameter :: pi=3.1415926535897932385d0
88real(8),
parameter :: g = 0.8098d0
89real(8),
parameter :: a = 2.804d0
90real(8),
parameter :: c = 0.7168d0
92real(8) s2,exlda, vxlda, eclda, vclda, X, Xs, Xss
93real(8) F, Fs, Fss, Hx, Hxs, Hxss, Hc, Hcs, Hcss
94real(8) zb, zbs, zbss, w
95real(8) n0b, n0bs, n0bss
96real(8) ln0b, ln0bs, ln0bss
97real(8) zbb, zbbc, zbbs, zbbss
98real(8) fxb, fxbs, fxbss
100if((rho <= 1.d-16))
then
116x = 1.0d0 - a*s2/(1.0d0 + a*s2)
119zb = (3.0d0/2.0d0*w)**(2.0d0/3.0d0)
120n0b = w/(2.0d0*pi**2*s**3)
121ln0b = -3.0d0/(2.0d0*pi)*(3.0d0*pi**2*n0b)**(1.0d0/3.0d0)
122zbbc = ((4.0d0/3.0d0)**(1.0d0/3.0d0)*2.0d0*pi/3.0d0)**4
123zbb = (zbbc*zb**2 + zb**4)**(1.0d0/4.0d0)
124fxb = -1.0d0/(ln0b*2.0d0*zbb)
125f = (c*s2 + 1.0d0)/(c*s2/fxb + 1.0d0)
127hx = x + (1.0d0 - x)*f
134hc = x + g*(1.0d0 - x)
142xs = -2.0d0*a*s/(1.0d0 + a*s2)**2
143xss = 2.0d0*a*(3.0d0*a*s2-1.0d0)/(1.0d0+a*s2)**3
146zbss = - zb*w*(5.0d0+2.0d0*w)/(2.0d0*s2*(1.0d0+w)**3)
147n0bs = sqrt(zb)*(-2.0d0*zb+s*zbs)/(2.0d0*pi**2*s2**2)
148n0bss = (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))
150ln0bs = -(3.0d0/pi)**(1.0d0/3.0d0)*n0bs/ &
151 (2.0d0*n0b**(2.0d0/3.0d0))
152ln0bss = (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))
154zbbs = zb*(zbbc+2*zb**2)*zbs/ &
155 (2.0d0*(zb**2*(zbbc+zb**2))**(3.0d0/4.0d0))
156zbbss = 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))
159fxbs = (zbb*ln0bs+ln0b*zbbs)/(2.0d0*ln0b**2*zbb**2)
160fxbss = (-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)
163fs = (c*s*(2.0d0*(fxb-1.0d0)*fxb + s*(1.0d0+c*s2)*fxbs))/ &
165fss = (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
170hxs = - (x - 1.0d0)*fs - (f - 1.0d0)*xs
171hxss = - 2.0d0*fs*xs - (x - 1.0d0)*fss - (f - 1.0d0)*xss
174vx = 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))
185vc = 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))
238real(8),
intent(in) :: n
239real(8),
intent(out) :: ec
240real(8),
intent(out) :: vc
242real(8),
parameter :: pi=3.1415926535897932385d0
243real(8),
parameter :: a01 = 0.21370d0
244real(8),
parameter :: b01 = 7.5957d0
245real(8),
parameter :: b02 = 3.5876d0
246real(8),
parameter :: b03 = 1.6382d0
247real(8),
parameter :: b04 = 0.49294d0
251real(8),
parameter :: A0 = 0.0310907d0
254real(8) Q0, Q1, Q1p, ecrs
255rsq = (3.0d0/(4.0d0*pi*n))**(1.0d0/6.0d0)
256ec = -2.0d0*a0*(1.0d0 + a01*rsq**2)* &
258 (2.0d0*a0*rsq*(b01 + rsq*(b02 + rsq*(b03 + b04*rsq)))))
259q0 = -2.0d0*a0*(1.0d0 + a01*rsq**2)
260q1 = 2.0d0*a0*rsq*(b01 + rsq*(b02 + rsq*(b03 + b04*rsq)))
261q1p = a0*(b01/rsq+2.0d0*b02+3.0d0*b03*rsq+4.0d0*b04*rsq**2)
262ecrs = -2.0d0*a0*a01*log(1.0d0 + 1.0d0/q1)-q0*q1p/(q1**2+q1)
263vc = ec - rsq**2/3.0d0*ecrs