The Elk Code
 
Loading...
Searching...
No Matches
potxcmt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2018 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine potxcmt(tsh,ias,xctype_,rhomt_,magmt_,taumt_,exmt_,ecmt_,vxcmt_, &
7 bxcmt_,wxcmt_)
8use modmain
9use modxcifc
10implicit none
11! arguments
12logical, intent(in) :: tsh
13integer, intent(in) :: ias,xctype_(3)
14real(8), intent(in) :: rhomt_(npmtmax,natmtot),magmt_(npmtmax,natmtot,ndmag)
15real(8), intent(in) :: taumt_(npmtmax,natmtot,nspinor)
16real(8), intent(out) :: exmt_(npmtmax,natmtot),ecmt_(npmtmax,natmtot)
17real(8), intent(out) :: vxcmt_(npmtmax,natmtot),bxcmt_(npmtmax,natmtot,ndmag)
18real(8), intent(out) :: wxcmt_(npmtmax,natmtot)
19! local variables
20integer ispn,idm,is
21integer nr,nri,n,i
22real(8) t0,t1,t2,t3,t4
23! allocatable arrays
24real(8), allocatable :: rho(:),rhoup(:),rhodn(:)
25real(8), allocatable :: gvrho(:,:),gvup(:,:),gvdn(:,:)
26real(8), allocatable :: grho(:),gup(:),gdn(:)
27real(8), allocatable :: g2rho(:),g2up(:),g2dn(:)
28real(8), allocatable :: g3rho(:),g3up(:),g3dn(:)
29real(8), allocatable :: grho2(:),gup2(:),gdn2(:),gupdn(:)
30real(8), allocatable :: ex(:),ec(:),vxc(:)
31real(8), allocatable :: vx(:),vxup(:),vxdn(:)
32real(8), allocatable :: vc(:),vcup(:),vcdn(:)
33real(8), allocatable :: mag(:,:),bxc(:,:),tau(:,:)
34real(8), allocatable :: dxdgr2(:),dxdgu2(:),dxdgd2(:),dxdgud(:)
35real(8), allocatable :: dcdgr2(:),dcdgu2(:),dcdgd2(:),dcdgud(:)
36real(8), allocatable :: dxdg2r(:),dxdg2u(:),dxdg2d(:)
37real(8), allocatable :: dcdg2r(:),dcdg2u(:),dcdg2d(:)
38real(8), allocatable :: dtdr(:),dtdru(:),dtdrd(:)
39real(8), allocatable :: wx(:),wxup(:),wxdn(:)
40real(8), allocatable :: wc(:),wcup(:),wcdn(:)
41is=idxis(ias)
42n=npmt(is)
43! allocate local arrays
44allocate(rho(n),ex(n),ec(n),vxc(n))
45if (any(xcgrad == [3,4,5,6])) allocate(tau(n,nspinor))
46if (spinpol) then
47 allocate(mag(n,3),bxc(n,3))
48end if
49if (spinpol) then
50 allocate(rhoup(n),rhodn(n))
51 allocate(vxup(n),vxdn(n),vcup(n),vcdn(n))
52 if (xcgrad == 1) then
53 allocate(grho(n),gup(n),gdn(n))
54 allocate(g2up(n),g2dn(n))
55 allocate(g3rho(n),g3up(n),g3dn(n))
56 else if (xcgrad == 2) then
57 allocate(g2up(n),g2dn(n))
58 allocate(gvup(n,3),gvdn(n,3))
59 allocate(gup2(n),gdn2(n),gupdn(n))
60 allocate(dxdgu2(n),dxdgd2(n),dxdgud(n))
61 allocate(dcdgu2(n),dcdgd2(n),dcdgud(n))
62 else if (any(xcgrad == [3,4,5,6])) then
63 allocate(g2up(n),g2dn(n))
64 allocate(gvup(n,3),gvdn(n,3))
65 allocate(gup2(n),gdn2(n),gupdn(n))
66 allocate(dxdgu2(n),dxdgd2(n),dxdgud(n))
67 allocate(dcdgu2(n),dcdgd2(n),dcdgud(n))
68 allocate(dxdg2u(n),dxdg2d(n))
69 allocate(dcdg2u(n),dcdg2d(n))
70 allocate(dtdru(n),dtdrd(n))
71 allocate(wxup(n),wxdn(n),wcup(n),wcdn(n))
72 end if
73else
74 allocate(vx(n),vc(n))
75 if (xcgrad == 1) then
76 allocate(grho(n),g2rho(n),g3rho(n))
77 else if (xcgrad == 2) then
78 allocate(g2rho(n),gvrho(n,3),grho2(n))
79 allocate(dxdgr2(n),dcdgr2(n))
80 else if (any(xcgrad == [3,4,5,6])) then
81 allocate(g2rho(n),gvrho(n,3),grho2(n))
82 allocate(dxdgr2(n),dcdgr2(n))
83 allocate(dxdg2r(n),dcdg2r(n))
84 allocate(dtdr(n),wx(n),wc(n))
85 end if
86end if
87nr=nrmt(is)
88nri=nrmti(is)
89if (tsh) then
90! convert the density to spherical coordinates
91 call rbsht(nr,nri,rhomt_(:,ias),rho)
92else
93 rho(1:n)=rhomt_(1:n,ias)
94end if
95! convert tau to spherical coordinates if required
96if (any(xcgrad == [3,4,5,6])) then
97 do ispn=1,nspinor
98 if (tsh) then
99 call rbsht(nr,nri,taumt_(:,ias,ispn),tau(:,ispn))
100 else
101 tau(1:n,ispn)=taumt_(1:n,ias,ispn)
102 end if
103 end do
104end if
105if (spinpol) then
106!------------------------!
107! spin-polarised !
108!------------------------!
109! magnetisation in spherical coordinates
110 do idm=1,ndmag
111 if (tsh) then
112 call rbsht(nr,nri,magmt_(:,ias,idm),mag(:,idm))
113 else
114 mag(1:n,idm)=magmt_(1:n,ias,idm)
115 end if
116 end do
117! use scaled spin exchange-correlation if required
118 if (tssxc) mag(1:n,1:ndmag)=mag(1:n,1:ndmag)*sxcscf
119 if (ncmag) then
120! non-collinear (use Kubler's trick)
121 if (xcgrad == 0) then
122! LSDA
123 do i=1,n
124! compute rhoup=(rho+|m|)/2 and rhodn=(rho-|m|)/2
125 t0=rho(i)
126 t1=sqrt(mag(i,1)**2+mag(i,2)**2+mag(i,3)**2)
127 rhoup(i)=0.5d0*(t0+t1)
128 rhodn(i)=0.5d0*(t0-t1)
129 end do
130 else
131! functionals which require gradients
132 do i=1,n
133 t0=rho(i)
134 t1=sqrt(mag(i,1)**2+mag(i,2)**2+mag(i,3)**2+dncgga)
135 rhoup(i)=0.5d0*(t0+t1)
136 rhodn(i)=0.5d0*(t0-t1)
137 end do
138 end if
139 else
140! collinear
141 do i=1,n
142! compute rhoup=(rho+m_z)/2 and rhodn=(rho-m_z)/2
143 t0=rho(i)
144 t1=mag(i,1)
145 rhoup(i)=0.5d0*(t0+t1)
146 rhodn(i)=0.5d0*(t0-t1)
147 end do
148 end if
149! call the exchange-correlation interface routine
150 if (xcgrad < 1) then
151! LSDA
152 call xcifc(xctype_,n,tempa=swidth,rhoup=rhoup,rhodn=rhodn,ex=ex,ec=ec, &
153 vxup=vxup,vxdn=vxdn,vcup=vcup,vcdn=vcdn)
154 else if (xcgrad == 1) then
155! GGA
156 call ggamt_sp_1(is,n,rhoup,rhodn,grho,gup,gdn,g2up,g2dn,g3rho,g3up,g3dn)
157 call xcifc(xctype_,n,rhoup=rhoup,rhodn=rhodn,grho=grho,gup=gup,gdn=gdn, &
158 g2up=g2up,g2dn=g2dn,g3rho=g3rho,g3up=g3up,g3dn=g3dn,ex=ex,ec=ec,vxup=vxup,&
159 vxdn=vxdn,vcup=vcup,vcdn=vcdn)
160 else if (xcgrad == 2) then
161! GGA
162 call ggamt_sp_2a(is,n,rhoup,rhodn,g2up,g2dn,gvup,gvdn,gup2,gdn2,gupdn)
163 call xcifc(xctype_,n,rhoup=rhoup,rhodn=rhodn,gup2=gup2,gdn2=gdn2, &
164 gupdn=gupdn,ex=ex,ec=ec,vxup=vxup,vxdn=vxdn,vcup=vcup,vcdn=vcdn, &
165 dxdgu2=dxdgu2,dxdgd2=dxdgd2,dxdgud=dxdgud,dcdgu2=dcdgu2,dcdgd2=dcdgd2, &
166 dcdgud=dcdgud)
167 call ggamt_sp_2b(is,n,g2up,g2dn,gvup,gvdn,vxup,vxdn,vcup,vcdn,dxdgu2, &
168 dxdgd2,dxdgud,dcdgu2,dcdgd2,dcdgud)
169 else if (any(xcgrad == [3,4,5])) then
170! meta-GGA with energy
171 call ggamt_sp_2a(is,n,rhoup,rhodn,g2up,g2dn,gvup,gvdn,gup2,gdn2,gupdn)
172! enforce the von Weizsacker lower bound
173 call k_vwlb(n,rhoup,gup2,tau(:,1))
174 call k_vwlb(n,rhodn,gdn2,tau(:,2))
175 call xcifc(xctype_,n,rhoup=rhoup,rhodn=rhodn,g2up=g2up,g2dn=g2dn,gup2=gup2,&
176 gdn2=gdn2,gupdn=gupdn,tauup=tau(:,1),taudn=tau(:,2),ex=ex,ec=ec,vxup=vxup,&
177 vxdn=vxdn,vcup=vcup,vcdn=vcdn,dxdgu2=dxdgu2,dxdgd2=dxdgd2,dxdgud=dxdgud, &
178 dcdgu2=dcdgu2,dcdgd2=dcdgd2,dcdgud=dcdgud,dxdg2u=dxdg2u,dxdg2d=dxdg2d, &
179 dcdg2u=dcdg2u,dcdg2d=dcdg2d,wxup=wxup,wxdn=wxdn,wcup=wcup,wcdn=wcdn)
180 call ggamt_sp_2b(is,n,g2up,g2dn,gvup,gvdn,vxup,vxdn,vcup,vcdn,dxdgu2, &
181 dxdgd2,dxdgud,dcdgu2,dcdgd2,dcdgud)
182! determine δτ(r')/δρ(r) using an approximate kinetic energy functional
183 if (any(xcgrad == [4,5])) then
184 call xcifc(ktype,n,rhoup=rhoup,rhodn=rhodn,g2up=g2up,g2dn=g2dn,gup2=gup2,&
185 gdn2=gdn2,tauup=tau(:,1),taudn=tau(:,2),dtdru=dtdru,dtdrd=dtdrd, &
186 dtdgu2=dxdgu2,dtdgd2=dxdgd2,dtdg2u=dxdg2u,dtdg2d=dxdg2d,wxup=dcdgu2, &
187 wxdn=dcdgd2)
188 call ggamt_4(is,n,gvup,vxup,vcup,wxup,wcup,dtdru,dxdgu2)
189 call ggamt_4(is,n,gvdn,vxdn,vcdn,wxdn,wcdn,dtdrd,dxdgd2)
190 if (kgrad == 3) then
191 call ggamt_3(is,n,vxup,vcup,wxup,wcup,dxdg2u)
192 call ggamt_3(is,n,vxdn,vcdn,wxdn,wcdn,dxdg2d)
193 end if
194 wxcmt_(1:n,ias)=0.5d0*(wxup(1:n)+wxdn(1:n)+wcup(1:n)+wcdn(1:n))
195 if (tsh) call rfshtip(nr,nri,wxcmt_(:,ias))
196 end if
197 else if (xcgrad == 6) then
198! potential-only meta-GGA
199 call ggamt_sp_2a(is,n,rhoup,rhodn,g2up,g2dn,gvup,gvdn,gup2,gdn2,gupdn)
200 call xcifc(xctype_,n,c_tb09=c_tb09,rhoup=rhoup,rhodn=rhodn,g2up=g2up, &
201 g2dn=g2dn,gup2=gup2,gdn2=gdn2,gupdn=gupdn,tauup=tau(:,1),taudn=tau(:,2), &
202 vxup=vxup,vxdn=vxdn,vcup=vcup,vcdn=vcdn,dxdgu2=dxdgu2,dxdgd2=dxdgd2, &
203 dxdgud=dxdgud,dcdgu2=dcdgu2,dcdgd2=dcdgd2,dcdgud=dcdgud,dxdg2u=dxdg2u, &
204 dxdg2d=dxdg2d,dcdg2u=dcdg2u,dcdg2d=dcdg2d,wxup=wxup,wxdn=wxdn,wcup=wcup, &
205 wcdn=wcdn)
206 ex(1:n)=0.d0; ec(1:n)=0.d0
207 wxcmt_(1:n,ias)=0.d0
208 end if
209! hybrid functionals
210 if (hybrid) then
211 t1=1.d0-hybridc
212! scale exchange part of energy
213 ex(1:n)=t1*ex(1:n)
214! scale exchange part of potential
215 vxup(1:n)=t1*vxup(1:n)
216 vxdn(1:n)=t1*vxdn(1:n)
217 end if
218 if (ncmag) then
219! non-collinear: locally spin rotate the exchange-correlation potential
220 do i=1,n
221 t1=vxup(i)+vcup(i)
222 t2=vxdn(i)+vcdn(i)
223 vxc(i)=0.5d0*(t1+t2)
224! determine the exchange-correlation magnetic field
225 t3=0.5d0*(t1-t2)
226! |m| = rhoup - rhodn
227 t4=rhoup(i)-rhodn(i)
228 if (abs(t4) > 1.d-8) t4=t3/t4
229 bxc(i,1:3)=mag(i,1:3)*t4
230 end do
231 else
232! collinear
233 do i=1,n
234 t1=vxup(i)+vcup(i)
235 t2=vxdn(i)+vcdn(i)
236 vxc(i)=0.5d0*(t1+t2)
237 bxc(i,1)=0.5d0*(t1-t2)
238 end do
239 end if
240! scale B_xc if required
241 if (tssxc) bxc(1:n,1:ndmag)=bxc(1:n,1:ndmag)*sxcscf
242 do idm=1,ndmag
243 if (tsh) then
244! convert field to spherical harmonics
245 call rfsht(nr,nri,bxc(:,idm),bxcmt_(:,ias,idm))
246 else
247 bxcmt_(1:n,ias,idm)=bxc(1:n,idm)
248 end if
249 end do
250else
251!--------------------------!
252! spin-unpolarised !
253!--------------------------!
254 if (xcgrad < 1) then
255 call xcifc(xctype_,n,tempa=swidth,rho=rho,ex=ex,ec=ec,vx=vx,vc=vc)
256 else if (xcgrad == 1) then
257 call ggamt_1(tsh,is,n,rhomt_(:,ias),grho,g2rho,g3rho)
258 call xcifc(xctype_,n,rho=rho,grho=grho,g2rho=g2rho,g3rho=g3rho,ex=ex,ec=ec,&
259 vx=vx,vc=vc)
260 else if (xcgrad == 2) then
261 call ggamt_2a(tsh,is,n,rhomt_(:,ias),g2rho,gvrho,grho2)
262 call xcifc(xctype_,n,rho=rho,grho2=grho2,ex=ex,ec=ec,vx=vx,vc=vc, &
263 dxdgr2=dxdgr2,dcdgr2=dcdgr2)
264 call ggamt_2b(is,n,g2rho,gvrho,vx,vc,dxdgr2,dcdgr2)
265 else if (any(xcgrad == [3,4,5])) then
266 call ggamt_2a(tsh,is,n,rhomt_(:,ias),g2rho,gvrho,grho2)
267! enforce the von Weizsacker lower bound
268 call k_vwlb(n,rho,grho2,tau)
269 call xcifc(xctype_,n,rho=rho,g2rho=g2rho,grho2=grho2,tau=tau,ex=ex,ec=ec, &
270 vx=vx,vc=vc,dxdgr2=dxdgr2,dcdgr2=dcdgr2,dxdg2r=dxdg2r,dcdg2r=dcdg2r,wx=wx,&
271 wc=wc)
272 call ggamt_2b(is,n,g2rho,gvrho,vx,vc,dxdgr2,dcdgr2)
273! determine δτ(r')/δρ(r) using an approximate kinetic energy functional
274 if (any(xcgrad == [4,5])) then
275 call xcifc(ktype,n,rho=rho,g2rho=g2rho,grho2=grho2,tau=tau,dtdr=dtdr, &
276 dtdgr2=dxdgr2,dtdg2r=dxdg2r,wx=dcdgr2)
277 call ggamt_4(is,n,gvrho,vx,vc,wx,wc,dtdr,dxdgr2)
278 if (kgrad == 3) then
279 call ggamt_3(is,n,vx,vc,wx,wc,dxdg2r)
280 end if
281 wxcmt_(1:n,ias)=wx(1:n)+wc(1:n)
282 if (tsh) call rfshtip(nr,nri,wxcmt_(:,ias))
283 end if
284 else if (xcgrad == 6) then
285 call ggamt_2a(tsh,is,n,rhomt_(:,ias),g2rho,gvrho,grho2)
286 call xcifc(xctype_,n,c_tb09=c_tb09,rho=rho,g2rho=g2rho,grho2=grho2,tau=tau,&
287 vx=vx,vc=vc,dxdgr2=dxdgr2,dcdgr2=dcdgr2,dxdg2r=dxdg2r,dcdg2r=dcdg2r, &
288 wx=wx,wc=wc)
289 ex(1:n)=0.d0; ec(1:n)=0.d0
290 wxcmt_(1:n,ias)=0.d0
291 end if
292! hybrid functionals
293 if (hybrid) then
294 t1=1.d0-hybridc
295! scale exchange part of energy
296 ex(1:n)=t1*ex(1:n)
297! scale exchange part of potential
298 vxc(1:n)=t1*vx(1:n)+vc(1:n)
299 else
300 vxc(1:n)=vx(1:n)+vc(1:n)
301 end if
302end if
303if (tsh) then
304! convert exchange and correlation energy densities to spherical harmonics
305 call rfsht(nr,nri,ex,exmt_(:,ias))
306 call rfsht(nr,nri,ec,ecmt_(:,ias))
307! convert exchange-correlation potential to spherical harmonics
308 call rfsht(nr,nri,vxc,vxcmt_(:,ias))
309else
310 exmt_(1:n,ias)=ex(1:n)
311 ecmt_(1:n,ias)=ec(1:n)
312 vxcmt_(1:n,ias)=vxc(1:n)
313end if
314deallocate(rho,ex,ec,vxc)
315if (any(xcgrad == [3,4,5,6])) deallocate(tau)
316if (spinpol) then
317 deallocate(mag,bxc)
318 deallocate(rhoup,rhodn,vxup,vxdn,vcup,vcdn)
319 if (xcgrad == 1) then
320 deallocate(grho,gup,gdn,g2up,g2dn,g3rho,g3up,g3dn)
321 else if (xcgrad == 2) then
322 deallocate(g2up,g2dn)
323 deallocate(gvup,gvdn)
324 deallocate(gup2,gdn2,gupdn)
325 deallocate(dxdgu2,dxdgd2,dxdgud)
326 deallocate(dcdgu2,dcdgd2,dcdgud)
327 else if (any(xcgrad == [3,4,5,6])) then
328 deallocate(g2up,g2dn)
329 deallocate(gvup,gvdn)
330 deallocate(gup2,gdn2,gupdn)
331 deallocate(dxdgu2,dxdgd2,dxdgud)
332 deallocate(dcdgu2,dcdgd2,dcdgud)
333 deallocate(dxdg2u,dxdg2d)
334 deallocate(dcdg2u,dcdg2d)
335 deallocate(dtdru,dtdrd)
336 deallocate(wxup,wxdn,wcup,wcdn)
337 end if
338else
339 deallocate(vx,vc)
340 if (xcgrad == 1) then
341 deallocate(grho,g2rho,g3rho)
342 else if (xcgrad == 2) then
343 deallocate(g2rho,gvrho,grho2)
344 deallocate(dxdgr2,dcdgr2)
345 else if (any(xcgrad == [3,4,5,6])) then
346 deallocate(g2rho,gvrho,grho2)
347 deallocate(dxdgr2,dcdgr2,dxdg2r,dcdg2r)
348 deallocate(dtdr,wx,wc)
349 end if
350end if
351end subroutine
352
subroutine ggamt_1(tsh, is, np, rho, grho, g2rho, g3rho)
Definition ggamt_1.f90:10
subroutine ggamt_2a(tsh, is, np, rho, g2rho, gvrho, grho2)
Definition ggamt_2a.f90:10
subroutine ggamt_2b(is, np, g2rho, gvrho, vx, vc, dxdgr2, dcdgr2)
Definition ggamt_2b.f90:10
subroutine ggamt_3(is, np, vx, vc, wx, wc, dtdg2r)
Definition ggamt_3.f90:7
subroutine ggamt_4(is, np, gvrho, vx, vc, wx, wc, dtdr, dtdgr2)
Definition ggamt_4.f90:7
subroutine ggamt_sp_1(is, np, rhoup, rhodn, grho, gup, gdn, g2up, g2dn, g3rho, g3up, g3dn)
subroutine ggamt_sp_2a(is, np, rhoup, rhodn, g2up, g2dn, gvup, gvdn, gup2, gdn2, gupdn)
subroutine ggamt_sp_2b(is, np, g2up, g2dn, gvup, gvdn, vxup, vxdn, vcup, vcdn, dxdgu2, dxdgd2, dxdgud, dcdgu2, dcdgd2, dcdgud)
pure subroutine k_vwlb(n, rho, grho2, tau)
Definition k_vwlb.f90:7
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8) hybridc
Definition modmain.f90:1151
logical ncmag
Definition modmain.f90:240
integer kgrad
Definition modmain.f90:610
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
logical spinpol
Definition modmain.f90:228
integer, dimension(3) ktype
Definition modmain.f90:606
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
logical hybrid
Definition modmain.f90:1149
real(8) swidth
Definition modmain.f90:892
real(8) c_tb09
Definition modmain.f90:678
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
logical tssxc
Definition modmain.f90:666
real(8) dncgga
Definition modmain.f90:604
integer xcgrad
Definition modmain.f90:602
real(8) sxcscf
Definition modmain.f90:668
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 potxcmt(tsh, ias, xctype_, rhomt_, magmt_, taumt_, exmt_, ecmt_, vxcmt_, bxcmt_, wxcmt_)
Definition potxcmt.f90:8
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition rbsht.f90:7
subroutine rfsht(nr, nri, rfmt1, rfmt2)
Definition rfsht.f90:7
subroutine rfshtip(nr, nri, rfmt)
Definition rfshtip.f90:7