The Elk Code
dpotxc.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2012 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 
6 subroutine dpotxc
7 use modmain
8 use modphonon
9 implicit none
10 ! local variables
11 integer idm,is,ias
12 integer nr,nri,nrc,nrci
13 integer ir,np
14 ! allocatable arrays
15 real(8), allocatable :: fxcmt(:,:,:,:),fxcir(:,:,:)
16 complex(8), allocatable :: dvmt(:),dbmt(:,:)
17 ! compute the exchange-correlation kernel
18 if (spinpol) then
19  allocate(fxcmt(npmtmax,natmtot,4,4),fxcir(ngtot,4,4))
20  call genspfxcr(.false.,fxcmt,fxcir)
21 else
22  allocate(fxcmt(npmtmax,natmtot,1,1),fxcir(ngtot,1,1))
23  call genfxcr(.false.,fxcmt,fxcir)
24 end if
25 allocate(dvmt(npmtmax))
26 if (spinpol) allocate(dbmt(npmtmax,3))
27 !---------------------------------------!
28 ! muffin-tin potential and field !
29 !---------------------------------------!
30 ! note: muffin-tin functions are in spherical coordinates
31 do ias=1,natmtot
32  is=idxis(ias)
33  nr=nrmt(is)
34  nri=nrmti(is)
35  nrc=nrcmt(is)
36  nrci=nrcmti(is)
37  np=npmt(is)
38 ! charge-charge contribution to potential derivative
39  dvmt(1:np)=fxcmt(1:np,ias,1,1)*drhomt(1:np,ias)
40 ! spin-polarised case
41  if (spinpol) then
42  if (ncmag) then
43 ! non-collinear
44 ! add charge-spin contribution to potential derivative
45  dvmt(1:np)=dvmt(1:np) &
46  +fxcmt(1:np,ias,1,2)*dmagmt(1:np,ias,1) &
47  +fxcmt(1:np,ias,1,3)*dmagmt(1:np,ias,2) &
48  +fxcmt(1:np,ias,1,4)*dmagmt(1:np,ias,3)
49 ! spin-charge contribution to B-field derivative
50  dbmt(1:np,1)=fxcmt(1:np,ias,1,2)*drhomt(1:np,ias)
51  dbmt(1:np,2)=fxcmt(1:np,ias,1,3)*drhomt(1:np,ias)
52  dbmt(1:np,3)=fxcmt(1:np,ias,1,4)*drhomt(1:np,ias)
53 ! add spin-spin contribution to B-field derivative
54 ! (note: f_xc is stored as an upper triangular matrix)
55  dbmt(1:np,1)=dbmt(1:np,1) &
56  +fxcmt(1:np,ias,2,2)*dmagmt(1:np,ias,1) &
57  +fxcmt(1:np,ias,2,3)*dmagmt(1:np,ias,2) &
58  +fxcmt(1:np,ias,2,4)*dmagmt(1:np,ias,3)
59  dbmt(1:np,2)=dbmt(1:np,2) &
60  +fxcmt(1:np,ias,2,3)*dmagmt(1:np,ias,1) &
61  +fxcmt(1:np,ias,3,3)*dmagmt(1:np,ias,2) &
62  +fxcmt(1:np,ias,3,4)*dmagmt(1:np,ias,3)
63  dbmt(1:np,3)=dbmt(1:np,3) &
64  +fxcmt(1:np,ias,2,4)*dmagmt(1:np,ias,1) &
65  +fxcmt(1:np,ias,3,4)*dmagmt(1:np,ias,2) &
66  +fxcmt(1:np,ias,4,4)*dmagmt(1:np,ias,3)
67  else
68 ! collinear
69 ! add charge-spin contribution to potential derivative
70  dvmt(1:np)=dvmt(1:np)+fxcmt(1:np,ias,1,4)*dmagmt(1:np,ias,1)
71 ! spin-charge contribution to B-field derivative
72  dbmt(1:np,1)=fxcmt(1:np,ias,1,4)*drhomt(1:np,ias)
73 ! add spin-spin contribution to B-field derivative
74  dbmt(1:np,1)=dbmt(1:np,1)+fxcmt(1:np,ias,4,4)*dmagmt(1:np,ias,1)
75  end if
76  end if
77 ! convert potential derivative to spherical harmonics
78  call zfsht(nr,nri,dvmt,dvsmt(:,ias))
79 ! convert magnetic field derivative to spherical harmonics on coarse mesh
80  do idm=1,ndmag
81  call zfmtftoc(nrc,nrci,dbmt(:,idm),dbsmt(:,ias,idm))
82  call zfshtip(nrc,nrci,dbsmt(:,ias,idm))
83  end do
84 end do
85 !------------------------------------------!
86 ! interstitial potential and field !
87 !------------------------------------------!
88 ! charge-charge contribution to potential derivative
89 do ir=1,ngtot
90  dvsir(ir)=fxcir(ir,1,1)*drhoir(ir)
91 end do
92 ! spin-polarised case
93 if (spinpol) then
94  if (ncmag) then
95 ! non-collinear
96  do ir=1,ngtot
97 ! add charge-spin contribution to potential derivative
98  dvsir(ir)=dvsir(ir) &
99  +fxcir(ir,1,2)*dmagir(ir,1) &
100  +fxcir(ir,1,3)*dmagir(ir,2) &
101  +fxcir(ir,1,4)*dmagir(ir,3)
102 ! spin-charge contribution to B-field derivative
103  dbsir(ir,1)=fxcir(ir,1,2)*drhoir(ir)
104  dbsir(ir,2)=fxcir(ir,1,3)*drhoir(ir)
105  dbsir(ir,3)=fxcir(ir,1,4)*drhoir(ir)
106 ! add spin-spin contribution to B-field derivative
107  dbsir(ir,1)=dbsir(ir,1) &
108  +fxcir(ir,2,2)*dmagir(ir,1) &
109  +fxcir(ir,2,3)*dmagir(ir,2) &
110  +fxcir(ir,2,4)*dmagir(ir,3)
111  dbsir(ir,2)=dbsir(ir,2) &
112  +fxcir(ir,2,3)*dmagir(ir,1) &
113  +fxcir(ir,3,3)*dmagir(ir,2) &
114  +fxcir(ir,3,4)*dmagir(ir,3)
115  dbsir(ir,3)=dbsir(ir,3) &
116  +fxcir(ir,2,4)*dmagir(ir,1) &
117  +fxcir(ir,3,4)*dmagir(ir,2) &
118  +fxcir(ir,4,4)*dmagir(ir,3)
119  end do
120  else
121 ! collinear
122  do ir=1,ngtot
123 ! add charge-spin contribution to potential derivative
124  dvsir(ir)=dvsir(ir)+fxcir(ir,1,4)*dmagir(ir,1)
125 ! spin-charge contribution to B-field derivative
126  dbsir(ir,1)=fxcir(ir,1,4)*drhoir(ir)
127 ! add spin-spin contribution to B-field derivative
128  dbsir(ir,1)=dbsir(ir,1)+fxcir(ir,4,4)*dmagir(ir,1)
129  end do
130  end if
131 end if
132 deallocate(fxcmt,fxcir,dvmt)
133 if (spinpol) deallocate(dbmt)
134 end subroutine
135 
subroutine genspfxcr(tsh, fxcmt, fxcir)
Definition: genspfxcr.f90:7
complex(8), dimension(:,:), pointer, contiguous dmagir
Definition: modphonon.f90:102
integer ngtot
Definition: modmain.f90:390
logical spinpol
Definition: modmain.f90:228
complex(8), dimension(:,:), pointer, contiguous drhomt
Definition: modphonon.f90:100
subroutine genfxcr(tsh, fxcmt, fxcir)
Definition: genfxcr.f90:7
complex(8), dimension(:), allocatable dvsir
Definition: modphonon.f90:112
integer ndmag
Definition: modmain.f90:238
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
complex(8), dimension(:,:,:), pointer, contiguous dbsmt
Definition: modphonon.f90:116
subroutine zfsht(nr, nri, zfmt1, zfmt2)
Definition: zfsht.f90:10
complex(8), dimension(:,:,:), pointer, contiguous dmagmt
Definition: modphonon.f90:102
pure subroutine zfmtftoc(nrc, nrci, zfmt, zfcmt)
Definition: zfmtftoc.f90:7
subroutine zfshtip(nr, nri, zfmt)
Definition: zfshtip.f90:7
complex(8), dimension(:,:), pointer, contiguous dvsmt
Definition: modphonon.f90:110
complex(8), dimension(:,:), pointer, contiguous dbsir
Definition: modphonon.f90:116
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
complex(8), dimension(:), pointer, contiguous drhoir
Definition: modphonon.f90:100
integer npmtmax
Definition: modmain.f90:216
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
logical ncmag
Definition: modmain.f90:240
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
subroutine dpotxc
Definition: dpotxc.f90:7
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150