The Elk Code
 
Loading...
Searching...
No Matches
zpotcoul.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!BOP
7! !ROUTINE: zpotcoul
8! !INTERFACE:
9subroutine zpotcoul(iash,nrmt_,nrmti_,npmt_,ld1,rl,ngridg_,igfft_,ngp,gpc, &
10 gclgp,ld2,jlgprmt,ylmgp,sfacgp,zrhoir,ld3,zvclmt,zvclir)
11! !USES:
12use modmain
13! !INPUT/OUTPUT PARAMETERS:
14! iash : if iash > 0 then the homogeous solution is stored for this
15! muffin-tin in the array zvclmt(:,natmtot+1) (in,integer)
16! nrmt_ : number of radial points for each species (in,integer(nspecies))
17! nrmti_ : number of radial points on inner part (in,integer(nspecies))
18! npmt_ : total number of points in muffin-tins (in,integer(nspecies))
19! ld1 : leading dimension (in,integer)
20! rl : r^l on radial mesh for each species
21! (in,real(ld1,-lmaxo-1:lmaxo+2,nspecies))
22! ngridg_ : G-vector grid sizes (in,integer(3))
23! igfft_ : map from G-vector index to FFT array (in,integer(*))
24! ngp : number of G+p-vectors (in,integer)
25! gpc : G+p-vector lengths (in,real(ngp))
26! gclgp : Coulomb Green's function in G+p-space (in,real(ngp))
27! ld2 : leading dimension (in,integer)
28! jlgprmt : spherical Bessel functions for evergy G+p-vector and muffin-tin
29! radius (in,real(0:lnpsd,ld2,nspecies))
30! ylmgp : spherical harmonics of the G+p-vectors (in,complex(lmmaxo,ngp))
31! sfacgp : structure factors of the G+p-vectors (in,complex(ld2,natmtot))
32! zrhoir : interstitial charge density (in,complex(*))
33! ld3 : leading dimension (in,integer)
34! zvclmt : muffin-tin Coulomb potential, with the contribution from the
35! isolated muffin-tin density precalculated and passed in
36! (inout,complex(ld3,*))
37! zvclir : interstitial Coulomb potential (out,complex(*))
38! !DESCRIPTION:
39! Calculates the Coulomb potential of a complex charge density by solving
40! Poisson's equation using the method of M. Weinert, {\it J. Math. Phys.}
41! {\bf 22}, 2433 (1981). First, the multipole moments of the muffin-tin charge
42! are determined for the $j$th atom of the $i$th species by
43! $$ q_{ij;lm}^{\rm MT}=\int_0^{R_i}r^{l+2}\rho_{ij;lm}(r)dr+z_{ij}Y_{00}
44! \,\delta_{l,0}\;, $$
45! where $R_i$ is the muffin-tin radius and $z_{ij}$ is a point charge located
46! at the atom center (usually the nuclear charge, which should be taken as
47! {\bf negative}). Next, the multipole moments of the continuation of the
48! interstitial density, $\rho^{\rm I}$, into the muffin-tin are found with
49! $$ q_{ij;lm}^{\rm I}=4\pi i^l R_i^{l+3}\sum_{\bf G}\frac{j_{l+1}(GR_i)}
50! {GR_i}\rho^{\rm I}({\bf G})\exp(i{\bf G}\cdot{\bf r}_{ij})Y_{lm}^*
51! (\hat{\bf G}), $$
52! remembering that
53! $$ \lim_{x\rightarrow 0}\frac{j_{l+n}(x)}{x^n}=\frac{1}{(2n+1)!!}
54! \delta_{l,0} $$
55! should be used for the case ${\bf G}=0$. A pseudocharge is now constructed
56! which is equal to the real density in the interstitial region and whose
57! multipoles are the difference between the real and interstitial muffin-tin
58! multipoles. This pseudocharge density is smooth in the sense that it can be
59! expanded in terms of the finite set of ${\bf G}$-vectors. In each muffin-tin
60! the pseudocharge has the form
61! $$ \rho_{ij}^{\rm P}({\bf r})=\rho^{\rm I}({\bf r}-{\bf r}_{ij})+\sum_{lm}
62! \rho_{ij;lm}^{\rm P}\frac{1}{R_i^{l+3}}\left(\frac{r}{R_i}\right)^l\left(1-
63! \frac{r^2}{R_i^2}\right)^{N_i}Y_{lm}(\hat{\bf r}) $$
64! where
65! $$ \rho_{ij;lm}^{\rm P}=\frac{(2l+2N_i+3)!!}{2^N_iN_i!(2l+1)!!}\left(
66! q_{ij;lm}^{\rm MT}-q_{ij;lm}^{\rm I}\right) $$
67! and $N_i\approx\frac{1}{4}R_iG_{\rm max}$ is generally a good choice.
68! The pseudocharge in reciprocal space is given by
69! $$ \rho^{\rm P}({\bf G})=\rho^{\rm I}({\bf G})+\sum_{ij;lm}2^{N_i}N_i!
70! \frac{4\pi(-i)^l}{\Omega R_i^l}\frac{j_{l+N_i+1}(GR_i)}{(GR_i)^{N_i+1}}
71! \rho_{ij;lm}^{\rm P}\exp(-i{\bf G}\cdot{\bf r}_{ij})Y_{lm}(\hat{\bf G}) $$
72! which may be used for solving Poisson's equation directly
73! $$ V^{\rm P}({\bf G})=\begin{cases}
74! 4\pi\frac{\rho^{\rm P}({\bf G})}{G^2} & G>0 \\
75! 0 & G=0 \end{cases} $$
76! The usual Green's function approach is then employed to determine the
77! potential in the muffin-tin sphere due to charge in the sphere. In other
78! words
79! $$ V_{ij;lm}^{\rm MT}(r)=\frac{4\pi}{2l+1}\left(\frac{1}{r^{l+1}}\int_0^r
80! \rho_{ij;lm}^{\rm MT}(r'){r'}^{l+2}dr'+r^l\int_r^{R_i}\frac{
81! \rho_{ij;lm}^{\rm MT}(r')}{{r'}^{l-1}}dr'\right)+\frac{1}{Y_{00}}
82! \frac{z_{ij}}{r}\delta_{l,0} $$
83! where the last term is the monopole arising from the point charge. All that
84! remains is to add the homogenous solution of Poisson's equation,
85! $$ V_{ij}^{\rm H}({\bf r})=\sum_{lm}V_{ij;lm}^{\rm H}\left(\frac{r}
86! {R_i}\right)^lY_{lm}(\hat{\bf r}), $$
87! to the muffin-tin potential so that it is continuous at the muffin-tin
88! boundary. Therefore the coefficients, $\rho_{ij;lm}^{\rm H}$, are given by
89! $$ V_{ij;lm}^{\rm H}=4\pi i^l\sum_{\bf G}j_{l}(Gr)V^{\rm P}({\bf G})
90! \exp(i{\bf G}\cdot{\bf r}_{ij})Y_{lm}^*(\hat{\bf G})-V_{ij;lm}^{\rm MT}
91! (R_i). $$
92! Finally note that the ${\bf G}$-vectors passed to the routine can represent
93! vectors with a non-zero offset, ${\bf G}+{\bf p}$ say, which is required for
94! calculating Coulomb matrix elements.
95!
96! !REVISION HISTORY:
97! Created April 2003 (JKD)
98!EOP
99!BOC
100implicit none
101! arguments
102integer, intent(in) :: iash,nrmt_(nspecies),nrmti_(nspecies),npmt_(nspecies)
103integer, intent(in) :: ld1
104real(8), intent(in) :: rl(ld1,-lmaxo-1:lmaxo+2,nspecies)
105integer, intent(in) :: ngridg_(3),igfft_(*),ngp
106real(8), intent(in) :: gpc(ngp),gclgp(ngp)
107integer, intent(in) :: ld2
108real(8), intent(in) :: jlgprmt(0:lnpsd,ld2,nspecies)
109complex(8), intent(in) :: ylmgp(lmmaxo,ngp),sfacgp(ld2,natmtot)
110complex(8), intent(in) :: zrhoir(*)
111integer, intent(in) :: ld3
112complex(8), intent(inout) :: zvclmt(ld3,*)
113complex(8), intent(out) :: zvclir(*)
114! local variables
115integer is,ia,ias
116integer nr,nri,iro,np
117integer l,lm,lma,lmb
118integer ig,jg,i,i0,i1
119real(8) t1,t2,t3
120complex(8) z1,z2
121! automatic arrays
122real(8) rl2(0:lmaxo)
123complex(8) qlm(lmmaxo,natmtot),zlm(lmmaxo),zhmt(ld3)
124! external functions
125real(8), external :: factn2
126! compute the multipole moments from the muffin-tin potentials
127t1=1.d0/fourpi
128do ias=1,natmtot
129 is=idxis(ias)
130 i=npmt_(is)-lmmaxo
131 do l=0,lmaxo
132 t2=t1*dble(2*l+1)*rmtl(l+1,is)
133 lma=l**2+1; lmb=lma+2*l
134 qlm(lma:lmb,ias)=t2*zvclmt(i+lma:i+lmb,ias)
135 end do
136end do
137! Fourier transform density to G-space and store in zvclir
138call zcopy(ngridg_(1)*ngridg_(2)*ngridg_(3),zrhoir,1,zvclir,1)
139call zfftifc(3,ngridg_,-1,zvclir)
140! subtract the multipole moments of the interstitial charge density
141do is=1,nspecies
142 rl2(0:lmaxo)=rmtl(2:lmaxo+2,is)
143 t1=rl2(0)*fourpi*y00
144 do ia=1,natoms(is)
145 ias=idxas(ia,is)
146 zlm(1:lmmaxo)=0.d0
147 do ig=1,ngp
148 jg=igfft_(ig)
149 if (gpc(ig) > epslat) then
150 z1=zvclir(jg)*sfacgp(ig,ias)/gpc(ig)
151 zlm(1)=zlm(1)+jlgprmt(1,ig,is)*t1*z1
152 do l=1,lmaxo
153 lma=l**2+1; lmb=lma+2*l
154 z2=jlgprmt(l+1,ig,is)*rl2(l)*z1
155 zlm(lma:lmb)=zlm(lma:lmb)+z2*conjg(ylmgp(lma:lmb,ig))
156 end do
157 else
158 t2=(fourpi/3.d0)*rmtl(3,is)*y00
159 zlm(1)=zlm(1)+t2*zvclir(jg)
160 end if
161 end do
162 qlm(1:lmmaxo,ias)=qlm(1:lmmaxo,ias)-zlm(1:lmmaxo)
163 end do
164end do
165! find the smooth pseudocharge within the muffin-tin whose multipoles are the
166! difference between the real muffin-tin and interstitial multipoles
167t1=factn2(2*lnpsd+1)/omega
168do ias=1,natmtot
169 is=idxis(ias)
170 do l=0,lmaxo
171 t2=t1/(factn2(2*l+1)*rmtl(l,is))
172 lma=l**2+1; lmb=lma+2*l
173 zlm(lma:lmb)=t2*qlm(lma:lmb,ias)
174 end do
175! add the pseudocharge and real interstitial densities in G-space
176 do ig=1,ngp
177 jg=igfft_(ig)
178 if (gpc(ig) > epslat) then
179 t2=gpc(ig)*rmt(is)
180 t3=1.d0/t2**lnpsd
181 z1=t3*fourpi*y00*zlm(1)
182 do l=1,lmaxo
183 lma=l**2+1; lmb=lma+2*l
184 t3=t3*t2
185 z1=z1+t3*sum(zlm(lma:lmb)*ylmgp(lma:lmb,ig))
186 end do
187 z2=jlgprmt(lnpsd,ig,is)*conjg(sfacgp(ig,ias))
188 zvclir(jg)=zvclir(jg)+z1*z2
189 else
190 t2=fourpi*y00/factn2(2*lnpsd+1)
191 zvclir(jg)=zvclir(jg)+t2*zlm(1)
192 end if
193 end do
194end do
195! solve Poisson's equation in G+p-space for the pseudocharge
196do ig=1,ngp
197 jg=igfft_(ig)
198 zvclir(jg)=gclgp(ig)*zvclir(jg)
199end do
200! match potentials at muffin-tin boundary by adding homogeneous solution
201do ias=1,natmtot
202 is=idxis(ias)
203 nr=nrmt_(is)
204 nri=nrmti_(is)
205 iro=nri+1
206 np=npmt_(is)
207! find the spherical harmonic expansion of the interstitial potential at the
208! muffin-tin radius
209 zlm(1:lmmaxo)=0.d0
210 do ig=1,ngp
211 z1=zvclir(igfft_(ig))*sfacgp(ig,ias)
212 zlm(1)=zlm(1)+(jlgprmt(0,ig,is)*fourpi*y00)*z1
213 do l=1,lmaxo
214 lma=l**2+1; lmb=lma+2*l
215 z2=jlgprmt(l,ig,is)*z1
216 zlm(lma:lmb)=zlm(lma:lmb)+z2*conjg(ylmgp(lma:lmb,ig))
217 end do
218 end do
219! calculate the homogenous solution
220 i=np-lmmaxo
221 do l=0,lmaxi
222 t1=1.d0/rmtl(l,is)
223 do lm=l**2+1,(l+1)**2
224 z1=t1*(zlm(lm)-zvclmt(i+lm,ias))
225 i1=lmmaxi*(nri-1)+lm
226 zhmt(lm:i1:lmmaxi)=z1*rl(1:nri,l,is)
227 i0=i1+lmmaxi
228 i1=lmmaxo*(nr-iro)+i0
229 zhmt(i0:i1:lmmaxo)=z1*rl(iro:nr,l,is)
230 end do
231 end do
232 do l=lmaxi+1,lmaxo
233 t1=1.d0/rmtl(l,is)
234 do lm=l**2+1,(l+1)**2
235 z1=t1*(zlm(lm)-zvclmt(i+lm,ias))
236 i0=lmmaxi*nri+lm
237 i1=lmmaxo*(nr-iro)+i0
238 zhmt(i0:i1:lmmaxo)=z1*rl(iro:nr,l,is)
239 end do
240 end do
241 zvclmt(1:np,ias)=zvclmt(1:np,ias)+zhmt(1:np)
242! store the homogenous solution if required
243 if (iash > 0) then
244 if (ias == iash) zvclmt(1:np,natmtot+1)=zhmt(1:np)
245 end if
246end do
247! Fourier transform interstitial potential to real-space
248call zfftifc(3,ngridg_,1,zvclir)
249end subroutine
250!EOC
251
elemental real(8) function factn2(n)
Definition factn2.f90:7
real(8), parameter y00
Definition modmain.f90:1233
real(8), dimension(:,:), allocatable rmtl
Definition modmain.f90:167
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(maxspecies) rmt
Definition modmain.f90:162
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer lmmaxi
Definition modmain.f90:207
real(8) omega
Definition modmain.f90:20
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), parameter fourpi
Definition modmain.f90:1231
real(8) epslat
Definition modmain.f90:24
integer lmaxi
Definition modmain.f90:205
subroutine zfftifc(nd, n, sgn, z)
subroutine zpotcoul(iash, nrmt_, nrmti_, npmt_, ld1, rl, ngridg_, igfft_, ngp, gpc, gclgp, ld2, jlgprmt, ylmgp, sfacgp, zrhoir, ld3, zvclmt, zvclir)
Definition zpotcoul.f90:11