The Elk Code
modxcifc.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 General Public License.
4 ! See the file COPYING for license details.
5 
6 module modxcifc
7 
8 use libxcifc
9 
10 contains
11 
12 !BOP
13 ! !ROUTINE: xcifc
14 ! !INTERFACE:
15 subroutine xcifc(xctype,n,c_tb09,tempa,rho,rhoup,rhodn,grho,gup,gdn,g2rho,g2up,&
16  g2dn,g3rho,g3up,g3dn,grho2,gup2,gdn2,gupdn,tau,tauup,taudn,ex,ec,vx,vc,vxup, &
17  vxdn,vcup,vcdn,dxdgr2,dxdgu2,dxdgd2,dxdgud,dcdgr2,dcdgu2,dcdgd2,dcdgud,dxdg2r,&
18  dxdg2u,dxdg2d,dcdg2r,dcdg2u,dcdg2d,wx,wxup,wxdn,wc,wcup,wcdn,dtdr,dtdru,dtdrd,&
19  dtdgr2,dtdgu2,dtdgd2,dtdg2r,dtdg2u,dtdg2d)
20 ! !INPUT/OUTPUT PARAMETERS:
21 ! xctype : type of exchange-correlation functional (in,integer(3))
22 ! n : number of density points (in,integer)
23 ! c_tb09 : Tran-Blaha '09 constant c (in,real,optional)
24 ! tempa : temperature in atomic units (in,real,optional)
25 ! rho : spin-unpolarised charge density (in,real(n),optional)
26 ! rhoup : spin-up charge density (in,real(n),optional)
27 ! rhodn : spin-down charge density (in,real(n),optional)
28 ! grho : |grad rho| (in,real(n),optional)
29 ! gup : |grad rhoup| (in,real(n),optional)
30 ! gdn : |grad rhodn| (in,real(n),optional)
31 ! g2rho : grad^2 rho (in,real(n),optional)
32 ! g2up : grad^2 rhoup (in,real(n),optional)
33 ! g2dn : grad^2 rhodn (in,real(n),optional)
34 ! g3rho : (grad rho).(grad |grad rho|) (in,real(n),optional)
35 ! g3up : (grad rhoup).(grad |grad rhoup|) (in,real(n),optional)
36 ! g3dn : (grad rhodn).(grad |grad rhodn|) (in,real(n),optional)
37 ! grho2 : |grad rho|^2 (in,real(n),optional)
38 ! gup2 : |grad rhoup|^2 (in,real(n),optional)
39 ! gdn2 : |grad rhodn|^2 (in,real(n),optional)
40 ! gupdn : (grad rhoup).(grad rhodn) (in,real(n),optional)
41 ! tau : kinetic energy density (in,real(n),optional)
42 ! tauup : spin-up kinetic energy density (in,real(n),optional)
43 ! taudn : spin-down kinetic energy density (in,real(n),optional)
44 ! ex : exchange energy density (out,real(n),optional)
45 ! ec : correlation energy density (out,real(n),optional)
46 ! vx : spin-unpolarised exchange potential (out,real(n),optional)
47 ! vc : spin-unpolarised correlation potential (out,real(n),optional)
48 ! vxup : spin-up exchange potential (out,real(n),optional)
49 ! vxdn : spin-down exchange potential (out,real(n),optional)
50 ! vcup : spin-up correlation potential (out,real(n),optional)
51 ! vcdn : spin-down correlation potential (out,real(n),optional)
52 ! dxdgr2 : de_x/d(|grad rho|^2) (out,real(n),optional)
53 ! dxdgu2 : de_x/d(|grad rhoup|^2) (out,real(n),optional)
54 ! dxdgd2 : de_x/d(|grad rhodn|^2) (out,real(n),optional)
55 ! dxdgud : de_x/d((grad rhoup).(grad rhodn)) (out,real(n),optional)
56 ! dcdgr2 : de_c/d(|grad rho|^2) (out,real(n),optional)
57 ! dcdgu2 : de_c/d(|grad rhoup|^2) (out,real(n),optional)
58 ! dcdgd2 : de_c/d(|grad rhodn|^2) (out,real(n),optional)
59 ! dcdgud : de_c/d((grad rhoup).(grad rhodn)) (out,real(n),optional)
60 ! dxdg2r : de_x/d(grad^2 rho) (out,real(n),optional)
61 ! dxdg2u : de_x/d(grad^2 rhoup) (out,real(n),optional)
62 ! dxdg2d : de_x/d(grad^2 rhodn) (out,real(n),optional)
63 ! dcdg2r : de_c/d(grad^2 rho) (out,real(n),optional)
64 ! dcdg2u : de_c/d(grad^2 rhoup) (out,real(n),optional)
65 ! dcdg2d : de_c/d(grad^2 rhodn) (out,real(n),optional)
66 ! wx : de_x/dtau (out,real(n),optional)
67 ! wxup : de_x/dtauup (out,real(n),optional)
68 ! wxdn : de_x/dtaudn (out,real(n),optional)
69 ! wc : de_c/dtau (out,real(n),optional)
70 ! wcup : de_c/dtauup (out,real(n),optional)
71 ! wcdn : de_c/dtaudn (out,real(n),optional)
72 ! dtdr : dtau/drho (out,real(n),optional)
73 ! dtdru : dtauup/drhoup (out,real(n),optional)
74 ! dtdrd : dtaudn/drhodn (out,real(n),optional)
75 ! dtdgr2 : dtau/d|grad rho|^2 (out,real(n),optional)
76 ! dtdgu2 : dtauup/d(|grad rhoup|^2) (out,real(n),optional)
77 ! dtdgd2 : dtaudn/d(|grad rhodn|^2) (out,real(n),optional)
78 ! dtdg2r : dtau/d(grad^2 rho) (out,real(n),optional)
79 ! dtdg2u : dtauup/d(grad^2 rhoup) (out,real(n),optional)
80 ! dtdg2d : dtaudn/d(grad^2 rhodn) (out,real(n),optional)
81 ! !DESCRIPTION:
82 ! Interface to the exchange-correlation routines. In the most general case
83 ! (meta-GGA), the exchange-correlation energy is given by
84 ! $$ E_{xc}[\rho^{\uparrow},\rho^{\downarrow}]=\int d^3r\,
85 ! \rho({\bf r})\,\varepsilon_{xc}(\rho^{\uparrow},\rho^{\downarrow},
86 ! |\nabla\rho|,|\nabla\rho^{\uparrow}|,|\nabla\rho^{\downarrow}|,
87 ! \nabla^2\rho^{\uparrow},\nabla^2\rho^{\downarrow},\tau), $$
88 ! where $\rho({\bf r})=\rho^{\uparrow}({\bf r})+\rho^{\downarrow}({\bf r})$ is
89 ! the density;
90 ! $$ \tau({\bf r})\equiv\sum_{i\;{\rm occ}}\nabla\psi({\bf r})\cdot
91 ! \nabla\psi({\bf r}) $$
92 ! is twice the spin-contracted kinetic energy density; and $\varepsilon_{xc}$
93 ! is the exchange-correlation energy per electron.
94 !
95 ! !REVISION HISTORY:
96 ! Created October 2002 (JKD)
97 !EOP
98 !BOC
99 implicit none
100 ! mandatory arguments
101 integer, intent(in) :: xctype(3),n
102 ! optional arguments
103 real(8), optional, intent(in) :: c_tb09,tempa
104 real(8), optional, intent(in) :: rho(n),rhoup(n),rhodn(n)
105 real(8), optional, intent(in) :: grho(n),gup(n),gdn(n)
106 real(8), optional, intent(in) :: g2rho(n),g2up(n),g2dn(n)
107 real(8), optional, intent(in) :: g3rho(n),g3up(n),g3dn(n)
108 real(8), optional, intent(in) :: grho2(n),gup2(n),gdn2(n),gupdn(n)
109 real(8), optional, intent(in) :: tau(n),tauup(n),taudn(n)
110 real(8), optional, intent(out) :: ex(n),ec(n),vx(n),vc(n)
111 real(8), optional, intent(out) :: vxup(n),vxdn(n),vcup(n),vcdn(n)
112 real(8), optional, intent(out) :: dxdgr2(n),dxdgu2(n),dxdgd2(n),dxdgud(n)
113 real(8), optional, intent(out) :: dxdg2r(n),dxdg2u(n),dxdg2d(n)
114 real(8), optional, intent(out) :: wx(n),wxup(n),wxdn(n)
115 real(8), optional, intent(out) :: dcdgr2(n),dcdgu2(n),dcdgd2(n),dcdgud(n)
116 real(8), optional, intent(out) :: dcdg2r(n),dcdg2u(n),dcdg2d(n)
117 real(8), optional, intent(out) :: wc(n),wcup(n),wcdn(n)
118 real(8), optional, intent(out) :: dtdr(n),dtdru(n),dtdrd(n)
119 real(8), optional, intent(out) :: dtdgr2(n),dtdgu2(n),dtdgd2(n)
120 real(8), optional, intent(out) :: dtdg2r(n),dtdg2u(n),dtdg2d(n)
121 ! local variables
122 real(8) kappa,mu,beta
123 ! allocatable arrays
124 real(8), allocatable :: ra(:,:)
125 if (n < 1) then
126  write(*,*)
127  write(*,'("Error(xcifc): n < 1 : ",I8)') n
128  write(*,*)
129  stop
130 end if
131 select case(abs(xctype(1)))
132 case(1)
133 ! No density-derived exchange-correlation energy or potential
134  if (present(ex)) ex(1:n)=0.d0
135  if (present(ec)) ec(1:n)=0.d0
136  if (present(vx)) vx(1:n)=0.d0
137  if (present(vc)) vc(1:n)=0.d0
138  if (present(vxup)) vxup(1:n)=0.d0
139  if (present(vxdn)) vxdn(1:n)=0.d0
140  if (present(vcup)) vcup(1:n)=0.d0
141  if (present(vcdn)) vcdn(1:n)=0.d0
142 case(2)
143 ! Perdew-Zunger parameterisation of Ceperley-Alder electron gas
144 ! J. Perdew and A. Zunger, Phys. Rev. B 23, 5048 (1981)
145 ! D. M. Ceperly and B. J. Alder, Phys. Rev. Lett. 45, 566 (1980)
146  if (present(rho).and.present(ex).and.present(ec).and.present(vx) &
147  .and.present(vc)) then
148  call xc_pzca(n,rho,ex,ec,vx,vc)
149  else
150  goto 10
151  end if
152 case(3)
153 ! Perdew-Wang parameterisation of the spin-polarised Ceperley-Alder electron gas
154 ! J. Perdew and Y. Wang, Phys. Rev. B 45, 13244 (1992)
155 ! D. M. Ceperly and B. J. Alder, Phys. Rev. Lett. 45, 566 (1980)
156  if (present(rhoup).and.present(rhodn).and.present(ex).and.present(ec) &
157  .and.present(vxup).and.present(vxdn).and.present(vcup) &
158  .and.present(vcdn)) then
159 ! spin-polarised density
160  call xc_pwca(n,rhoup,rhodn,ex,ec,vxup,vxdn,vcup,vcdn)
161  else if (present(rho).and.present(ex).and.present(ec).and.present(vx) &
162  .and.present(vc)) then
163 ! divide spin-unpolarised density into up and down
164  allocate(ra(n,2))
165  ra(1:n,1)=0.5d0*rho(1:n)
166  call xc_pwca(n,ra,ra,ex,ec,vx,ra(:,2),vc,ra(:,2))
167  deallocate(ra)
168  else
169  goto 10
170  end if
171 case(4)
172 ! X-alpha approximation
173 ! J. C. Slater, Phys. Rev. 81, 385 (1951)
174  if (present(rho).and.present(ex).and.present(ec).and.present(vx) &
175  .and.present(vc)) then
176  call xc_xalpha(n,rho,ex,vx)
177 ! set correlation energy and potential to zero
178  ec(1:n)=0.d0
179  vc(1:n)=0.d0
180  else
181  goto 10
182  end if
183 case(5)
184 ! U. von Barth and L. Hedin parameterisation of LSDA
185 ! J. Phys. C, 5, 1629 (1972)
186  if (present(rhoup).and.present(rhodn).and.present(ex).and.present(ec) &
187  .and.present(vxup).and.present(vxdn).and.present(vcup) &
188  .and.present(vcdn)) then
189 ! spin-polarised density
190  call xc_vbh(n,rhoup,rhodn,ex,ec,vxup,vxdn,vcup,vcdn)
191  else if (present(rho).and.present(ex).and.present(ec).and.present(vx) &
192  .and.present(vc)) then
193 ! divide spin-unpolarised density into up and down
194  allocate(ra(n,2))
195  ra(1:n,1)=0.5d0*rho(1:n)
196  call xc_vbh(n,ra,ra,ex,ec,vx,ra(:,2),vc,ra(:,2))
197  deallocate(ra)
198  else
199  goto 10
200  end if
201 case(20,21,22)
202 ! original PBE kappa
203  kappa=0.804d0
204  if (xctype(1) == 21) then
205 ! Zhang-Yang kappa
206  kappa=1.245d0
207  end if
208 ! original PBE mu and beta
209  mu=0.2195149727645171d0
210  beta=0.06672455060314922d0
211  if (xctype(1) == 22) then
212 ! PBEsol parameters
213  mu=10.d0/81.d0
214  beta=0.046d0
215  end if
216 ! Perdew-Burke-Ernzerhof generalised gradient approximation
217 ! Phys. Rev. Lett. 77, 3865 (1996); 78, 1396(E) (1997)
218 ! Revised PBE, Zhang-Yang, Phys. Rev. Lett. 80, 890 (1998)
219  if (present(rhoup).and.present(rhodn).and.present(grho).and.present(gup) &
220  .and.present(gdn).and.present(g2up).and.present(g2dn).and.present(g3rho) &
221  .and.present(g3up).and.present(g3dn).and.present(ex).and.present(ec) &
222  .and.present(vxup).and.present(vxdn).and.present(vcup) &
223  .and.present(vcdn)) then
224  call xc_pbe(n,kappa,mu,beta,rhoup,rhodn,grho,gup,gdn,g2up,g2dn,g3rho,g3up, &
225  g3dn,ex,ec,vxup,vxdn,vcup,vcdn)
226  else if (present(rho).and.present(grho).and.present(g2rho) &
227  .and.present(g3rho).and.present(ex).and.present(ec).and.present(vx) &
228  .and.present(vc)) then
229  allocate(ra(n,5))
230  ra(1:n,1)=0.5d0*rho(1:n)
231  ra(1:n,2)=0.5d0*grho(1:n)
232  ra(1:n,3)=0.5d0*g2rho(1:n)
233  ra(1:n,4)=0.25d0*g3rho(1:n)
234  call xc_pbe(n,kappa,mu,beta,ra,ra,grho,ra(:,2),ra(:,2),ra(:,3),ra(:,3), &
235  g3rho,ra(:,4),ra(:,4),ex,ec,vx,ra(:,5),vc,ra(:,5))
236  deallocate(ra)
237  else
238  goto 10
239  end if
240 case(26)
241 ! Wu-Cohen exchange with PBE correlation generalised gradient functional
242 ! Zhigang Wu and R. E. Cohen, Phys. Rev. B 73, 235116 (2006)
243  if (present(rho).and.present(grho).and.present(g2rho).and.present(g3rho) &
244  .and.present(ex).and.present(ec).and.present(vx).and.present(vc)) then
245  call xc_wc06(n,rho,grho,g2rho,g3rho,ex,ec,vx,vc)
246  else
247  goto 10
248  end if
249 case(30)
250 ! Armiento-Mattsson generalised gradient functional
251 ! R. Armiento and A. E. Mattsson, Phys. Rev. B 72, 085108 (2005)
252  if (present(rho).and.present(grho).and.present(g2rho).and.present(g3rho) &
253  .and.present(ex).and.present(ec).and.present(vx).and.present(vc)) then
254  call xc_am05(n,rho,grho,g2rho,g3rho,ex,ec,vx,vc)
255  else
256  goto 10
257  end if
258 case(50)
259 ! Thomas-Fermi kinetic energy functional derivative δτ(r')/δρ(r)
260  if (present(rhoup).and.present(rhodn).and.present(dtdru) &
261  .and.present(dtdrd)) then
262  call k_tf_sp(n,rhoup,rhodn,dtdru,dtdrd)
263  if (present(dtdgu2)) dtdgu2(1:n)=0.d0
264  if (present(dtdgd2)) dtdgd2(1:n)=0.d0
265  else if (present(rho).and.present(dtdr)) then
266  call k_tf(n,rho,dtdr)
267  if (present(dtdgr2)) dtdgr2(1:n)=0.d0
268  else
269  goto 10
270  end if
271 case(51)
272 ! Thomas-Fermi kinetic energy functional derivative scaled by the ratio of the
273 ! exact τ to the TF τ
274  if (present(rhoup).and.present(rhodn).and.present(tauup).and.present(taudn) &
275  .and.present(dtdru).and.present(dtdrd)) then
276  call k_tfsc_sp(n,rhoup,rhodn,tauup,taudn,dtdru,dtdrd)
277  if (present(dtdgu2)) dtdgu2(1:n)=0.d0
278  if (present(dtdgd2)) dtdgd2(1:n)=0.d0
279  else if (present(rho).and.present(tau).and.present(dtdr)) then
280  call k_tfsc(n,rho,tau,dtdr)
281  if (present(dtdgr2)) dtdgr2(1:n)=0.d0
282  else
283  goto 10
284  end if
285 case(52)
286 ! Thomas-Fermi-von Weizsacker kinetic energy functional derivative
287  if (present(rhoup).and.present(rhodn).and.present(gup2).and.present(gdn2) &
288  .and.present(dtdru).and.present(dtdrd).and.present(dtdgu2) &
289  .and.present(dtdgd2)) then
290  call k_tfvw_sp(n,rhoup,rhodn,gup2,gdn2,dtdru,dtdrd,dtdgu2,dtdgd2)
291  else if (present(rho).and.present(grho2).and.present(dtdr) &
292  .and.present(dtdgr2)) then
293  call k_tfvw(n,rho,grho2,dtdr,dtdgr2)
294  else
295  goto 10
296  end if
297 case(100)
298 ! libxc library functionals
299  if (present(rhoup).and.present(rhodn).and.present(g2up).and.present(g2dn) &
300  .and.present(gup2).and.present(gdn2).and.present(gupdn).and.present(tauup) &
301  .and.present(taudn).and.present(vxup).and.present(vxdn).and.present(vcup) &
302  .and.present(vcdn).and.present(dxdgu2).and.present(dxdgd2) &
303  .and.present(dxdgud).and.present(dcdgu2).and.present(dcdgd2) &
304  .and.present(dcdgud).and.present(dxdg2u).and.present(dxdg2d) &
305  .and.present(dcdg2u).and.present(dcdg2d).and.present(wxup).and.present(wxdn)&
306  .and.present(wcup).and.present(wcdn)) then
307 ! spin-polarised meta-GGA
308  if (present(ex).and.present(ec)) then
309 ! energy functional
310  call xcifc_libxc(xctype,n,rhoup=rhoup,rhodn=rhodn,g2up=g2up,g2dn=g2dn, &
311  gup2=gup2,gdn2=gdn2,gupdn=gupdn,tauup=tauup,taudn=taudn,ex=ex,ec=ec, &
312  vxup=vxup,vxdn=vxdn,vcup=vcup,vcdn=vcdn,dxdgu2=dxdgu2,dxdgd2=dxdgd2, &
313  dxdgud=dxdgud,dcdgu2=dcdgu2,dcdgd2=dcdgd2,dcdgud=dcdgud,dxdg2u=dxdg2u, &
314  dxdg2d=dxdg2d,dcdg2u=dcdg2u,dcdg2d=dcdg2d,wxup=wxup,wxdn=wxdn,wcup=wcup,&
315  wcdn=wcdn)
316  else
317 ! potential-only functional
318  call xcifc_libxc(xctype,n,c_tb09=c_tb09,rhoup=rhoup,rhodn=rhodn, &
319  g2up=g2up,g2dn=g2dn,gup2=gup2,gdn2=gdn2,gupdn=gupdn,tauup=tauup, &
320  taudn=taudn,vxup=vxup,vxdn=vxdn,vcup=vcup,vcdn=vcdn,dxdgu2=dxdgu2, &
321  dxdgd2=dxdgd2,dxdgud=dxdgud,dcdgu2=dcdgu2,dcdgd2=dcdgd2,dcdgud=dcdgud, &
322  dxdg2u=dxdg2u,dxdg2d=dxdg2d,dcdg2u=dcdg2u,dcdg2d=dcdg2d,wxup=wxup, &
323  wxdn=wxdn,wcup=wcup,wcdn=wcdn)
324  end if
325  else if (present(rhoup).and.present(rhodn).and.present(gup2) &
326  .and.present(gdn2).and.present(gupdn).and.present(ex).and.present(ec) &
327  .and.present(vxup).and.present(vxdn).and.present(vcup).and.present(vcdn) &
328  .and.present(dxdgu2).and.present(dxdgd2).and.present(dxdgud) &
329  .and.present(dcdgu2).and.present(dcdgd2).and.present(dcdgud)) then
330 ! spin-polarised GGA
331  call xcifc_libxc(xctype,n,rhoup=rhoup,rhodn=rhodn,gup2=gup2,gdn2=gdn2, &
332  gupdn=gupdn,ex=ex,ec=ec,vxup=vxup,vxdn=vxdn,vcup=vcup,vcdn=vcdn, &
333  dxdgu2=dxdgu2,dxdgd2=dxdgd2,dxdgud=dxdgud,dcdgu2=dcdgu2,dcdgd2=dcdgd2, &
334  dcdgud=dcdgud)
335  else if (present(rhoup).and.present(rhodn).and.present(g2up) &
336  .and.present(g2dn).and.present(gup2).and.present(gdn2).and.present(tauup) &
337  .and.present(taudn).and.present(dtdru).and.present(dtdrd) &
338  .and.present(dtdg2u).and.present(dtdg2d).and.present(dtdgu2) &
339  .and.present(dtdgd2).and.present(wxup).and.present(wxdn)) then
340 ! spin-polarised kinetic energy functional
341  call xcifc_libxc(xctype,n,rhoup=rhoup,rhodn=rhodn,g2up=g2up,g2dn=g2dn, &
342  gup2=gup2,gdn2=gdn2,tauup=tauup,taudn=taudn,vxup=dtdru,vxdn=dtdrd, &
343  dxdg2u=dtdg2u,dxdg2d=dtdg2d,dxdgu2=dtdgu2,dxdgd2=dtdgd2,wxup=wxup, &
344  wxdn=wxdn)
345  else if (present(rhoup).and.present(rhodn).and.present(ex).and.present(ec) &
346  .and.present(vxup).and.present(vxdn).and.present(vcup) &
347  .and.present(vcdn)) then
348 ! spin-polarised LSDA
349  call xcifc_libxc(xctype,n,tempa=tempa,rhoup=rhoup,rhodn=rhodn,ex=ex,ec=ec, &
350  vxup=vxup,vxdn=vxdn,vcup=vcup,vcdn=vcdn)
351  else if (present(rho).and.present(g2rho).and.present(grho2).and.present(tau) &
352  .and.present(vx).and.present(vc).and.present(dxdgr2).and.present(dcdgr2) &
353  .and.present(dxdg2r).and.present(dcdg2r).and.present(wx).and.present(wc))then
354 ! spin-unpolarised meta-GGA
355  if (present(ex).and.present(ec)) then
356 ! energy functional
357  call xcifc_libxc(xctype,n,rho=rho,g2rho=g2rho,grho2=grho2,tau=tau,ex=ex, &
358  ec=ec,vx=vx,vc=vc,dxdgr2=dxdgr2,dcdgr2=dcdgr2,dxdg2r=dxdg2r, &
359  dcdg2r=dcdg2r,wx=wx,wc=wc)
360  else
361 ! potential-only functional
362  call xcifc_libxc(xctype,n,c_tb09=c_tb09,rho=rho,g2rho=g2rho,grho2=grho2, &
363  tau=tau,vx=vx,vc=vc,dxdgr2=dxdgr2,dcdgr2=dcdgr2,dxdg2r=dxdg2r, &
364  dcdg2r=dcdg2r,wx=wx,wc=wc)
365  end if
366  else if (present(rho).and.present(grho2).and.present(ex).and.present(ec) &
367  .and.present(vx).and.present(vc).and.present(dxdgr2) &
368  .and.present(dcdgr2)) then
369 ! spin-unpolarised GGA
370  call xcifc_libxc(xctype,n,rho=rho,grho2=grho2,ex=ex,ec=ec,vx=vx,vc=vc, &
371  dxdgr2=dxdgr2,dcdgr2=dcdgr2)
372  else if (present(rho).and.present(g2rho).and.present(grho2).and.present(tau) &
373  .and.present(dtdr).and.present(dtdgr2).and.present(dtdg2r) &
374  .and.present(wx)) then
375 ! spin-unpolarised kinetic energy functional
376  call xcifc_libxc(xctype,n,rho=rho,g2rho=g2rho,grho2=grho2,tau=tau,vx=dtdr, &
377  dxdgr2=dtdgr2,dxdg2r=dtdg2r,wx=wx)
378  else if (present(rho).and.present(ex).and.present(ec).and.present(vx) &
379  .and.present(vc)) then
380 ! spin-unpolarised LDA
381  call xcifc_libxc(xctype,n,tempa=tempa,rho=rho,ex=ex,ec=ec,vx=vx,vc=vc)
382  else
383  goto 10
384  end if
385 case default
386  write(*,*)
387  write(*,'("Error(xcifc): xctype not defined : ",I8)') xctype(1)
388  write(*,*)
389  stop
390 end select
391 ! set exchange potential to zero for EXX
392 if (xctype(1) <= -2) then
393  if (present(vx)) vx(1:n)=0.d0
394  if (present(vxup)) vxup(1:n)=0.d0
395  if (present(vxdn)) vxdn(1:n)=0.d0
396 end if
397 return
398 10 continue
399 write(*,*)
400 write(*,'("Error(xcifc): missing arguments for exchange-correlation type ",&
401  &3I6)') xctype(:)
402 write(*,*)
403 stop
404 end subroutine
405 !EOC
406 
407 !BOP
408 ! !ROUTINE: getxcdata
409 ! !INTERFACE:
410 subroutine getxcdata(xctype,xcdescr,xcspin,xcgrad,hybrid,hybridc)
411 ! !INPUT/OUTPUT PARAMETERS:
412 ! xctype : type of exchange-correlation functional (in,integer(3))
413 ! xcdescr : description of functional (out,character(264))
414 ! xcspin : spin treatment (out,integer)
415 ! xcgrad : gradient treatment (out,integer)
416 ! hybrid : .true. if functional a hybrid (out,logical)
417 ! hybridc : hybrid exact exchange mixing coefficient (out,real(8))
418 ! !DESCRIPTION:
419 ! Returns data on the exchange-correlation functional labeled by {\tt xctype}.
420 ! The character array {\tt xcdescr} contains a short description of the
421 ! functional including journal references. The variable {\tt xcspin} is set to
422 ! 1 or 0 for spin-polarised or -unpolarised functionals, respectively. For
423 ! functionals which require the gradients of the density {\tt xcgrad} is set
424 ! to 1, otherwise it is set to 0.
425 !
426 ! !REVISION HISTORY:
427 ! Created October 2002 (JKD)
428 !EOP
429 !BOC
430 implicit none
431 ! arguments
432 integer, intent(in) :: xctype(3)
433 character(264), intent(out) :: xcdescr
434 integer, intent(out) :: xcspin,xcgrad
435 logical, intent(out) :: hybrid
436 real(8), intent(out) :: hybridc
437 select case(abs(xctype(1)))
438 case(1)
439  xcdescr='No density-derived exchange-correlation energy or potential'
440 ! spin-polarisation or gradient status not required
441  xcspin=-1
442  xcgrad=-1
443 case(2)
444  xcdescr='Perdew-Zunger/Ceperley-Alder, Phys. Rev. B 23, 5048 (1981)'
445  xcspin=0
446  xcgrad=0
447 case(3)
448  xcdescr='Perdew-Wang/Ceperley-Alder, Phys. Rev. B 45, 13244 (1992)'
449  xcspin=1
450  xcgrad=0
451 case(4)
452  xcdescr='X-alpha approximation, J. C. Slater, Phys. Rev. 81, 385 (1951)'
453  xcspin=0
454  xcgrad=0
455 case(5)
456  xcdescr='von Barth-Hedin, J. Phys. C 5, 1629 (1972)'
457  xcspin=1
458  xcgrad=0
459 case(20)
460  xcdescr='Perdew-Burke-Ernzerhof, Phys. Rev. Lett. 77, 3865 (1996)'
461  xcspin=1
462  xcgrad=1
463 case(21)
464  xcdescr='Revised PBE, Zhang-Yang, Phys. Rev. Lett. 80, 890 (1998)'
465  xcspin=1
466  xcgrad=1
467 case(22)
468  xcdescr='PBEsol, Phys. Rev. Lett. 100, 136406 (2008)'
469  xcspin=1
470  xcgrad=1
471 case(26)
472  xcdescr='Wu-Cohen exchange + PBE correlation, Phys. Rev. B 73, 235116 (2006)'
473  xcspin=0
474  xcgrad=1
475 case(30)
476  xcdescr='Armiento-Mattsson functional, Phys. Rev. B 72, 85108 (2005)'
477  xcspin=0
478  xcgrad=1
479 case(50)
480  xcdescr='Thomas-Fermi kinetic energy functional'
481  xcspin=-1
482  xcgrad=0
483 case(51)
484  xcdescr='Thomas-Fermi kinetic energy functional scaled by the ratio of the &
485  &exact τ to the TF τ'
486  xcspin=-1
487  xcgrad=0
488 case(52)
489  xcdescr='Thomas-Fermi-von Weizsacker kinetic energy functional'
490  xcspin=-1
491  xcgrad=1
492 case(100)
493 ! libxc library functionals
494  call xcdata_libxc(xctype,xcdescr,xcspin,xcgrad,hybrid,hybridc)
495 case default
496  write(*,*)
497  write(*,'("Error(getxcdata): xctype not defined : ",I8)') xctype(1)
498  write(*,*)
499  stop
500 end select
501 end subroutine
502 !EOC
503 
504 end module
505 
subroutine k_tf_sp(n, rhoup, rhodn, dtdru, dtdrd)
Definition: k_tf_sp.f90:7
subroutine xc_pbe(n, kappa, mu, beta, rhoup, rhodn, grho, gup, gdn, g2up, g2dn, g3rho, g3up, g3dn, ex, ec, vxup, vxdn, vcup, vcdn)
Definition: xc_pbe.f90:9
subroutine xc_xalpha(n, rho, exc, vxc)
Definition: xc_xalpha.f90:10
subroutine k_tfsc(n, rho, tau, dtdr)
Definition: k_tfsc.f90:7
subroutine k_tf(n, rho, dtdr)
Definition: k_tf.f90:7
subroutine xcifc_libxc(xctype, n, c_tb09, tempa, rho, rhoup, rhodn, g2rho, g2up, g2dn, grho2, gup2, gdn2, gupdn, tau, tauup, taudn, ex, ec, vx, vc, vxup, vxdn, vcup, vcdn, dxdgr2, dxdgu2, dxdgd2, dxdgud, dcdgr2, dcdgu2, dcdgd2, dcdgud, dxdg2r, dxdg2u, dxdg2d, dcdg2r, dcdg2u, dcdg2d, wx, wxup, wxdn, wc, wcup, wcdn)
subroutine xcdata_libxc(xctype, xcdescr, xcspin, xcgrad, hybrid, hybridc)
subroutine getxcdata(xctype, xcdescr, xcspin, xcgrad, hybrid, hybridc)
Definition: modxcifc.f90:411
subroutine k_tfsc_sp(n, rhoup, rhodn, tauup, taudn, dtdru, dtdrd)
Definition: k_tfsc_sp.f90:7
subroutine k_tfvw_sp(n, rhoup, rhodn, gup2, gdn2, dtdru, dtdrd, dtdgu2, dtdgd2)
Definition: k_tfvw_sp.f90:10
subroutine xc_wc06(n, rho, grho, g2rho, g3rho, ex, ec, vx, vc)
Definition: xc_wc06.f90:7
subroutine xc_pwca(n, rhoup, rhodn, ex, ec, vxup, vxdn, vcup, vcdn)
Definition: xc_pwca.f90:10
subroutine xc_am05(n, rho, grho, g2rho, g3rho, ex, ec, vx, vc)
Definition: xc_am05.f90:10
subroutine xc_vbh(n, rhoup, rhodn, ex, ec, vxup, vxdn, vcup, vcdn)
Definition: xc_vbh.f90:10
subroutine k_tfvw(n, rho, grho2, dtdr, dtdgr2)
Definition: k_tfvw.f90:10
subroutine xcifc(xctype, n, c_tb09, tempa, rho, rhoup, rhodn, grho, gup, gdn, g2rho, g2up, g2dn, g3rho, g3up, g3dn, grho2, gup2, gdn2, gupdn, tau, tauup, taudn, ex, ec, vx, vc, vxup, vxdn, vcup, vcdn, dxdgr2, dxdgu2, dxdgd2, dxdgud, dcdgr2, dcdgu2, dcdgd2, dcdgud, dxdg2r, dxdg2u, dxdg2d, dcdg2r, dcdg2u, dcdg2d, wx, wxup, wxdn, wc, wcup, wcdn, dtdr, dtdru, dtdrd, dtdgr2, dtdgu2, dtdgd2, dtdg2r, dtdg2u, dtdg2d)
Definition: modxcifc.f90:20
subroutine xc_pzca(n, rho, ex, ec, vx, vc)
Definition: xc_pzca.f90:10