The Elk Code
 
Loading...
Searching...
No Matches
gradzfmt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2009 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: gradzfmt
8! !INTERFACE:
9subroutine gradzfmt(nr,nri,ri,wcr,zfmt,ld,gzfmt)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! nr : number of radial mesh points (in,integer)
14! nri : number of points on inner part of muffin-tin (in,integer)
15! ri : 1/r on the radial mesh (in,real(nr))
16! wcr : weights for spline coefficients on radial mesh (in,real(12,nr))
17! zfmt : complex muffin-tin function (in,complex(*))
18! ld : leading dimension (in,integer)
19! gzfmt : gradient of zfmt (out,complex(ld,3))
20! !DESCRIPTION:
21! Calculates the gradient of a complex muffin-tin function. In other words,
22! given the spherical harmonic expansion coefficients, $f_{lm}(r)$, of a
23! function $f({\bf r})$, the routine returns ${\bf F}_{lm}$ where
24! $$ \sum_{lm}{\bf F}_{lm}(r)Y_{lm}(\hat{\bf r})=\nabla f({\bf r}). $$
25! This is done using the gradient formula (see, for example, V. Devanathan,
26! {\em Angular Momentum Techniques In Quantum Mechanics})
27! \begin{align*}
28! \nabla f_{lm}(r)Y_{lm}(\hat{\bf r})&=-\sqrt{\frac{l+1}{2l+1}}
29! \left(\frac{d}{dr}-\frac{l}{r}\right)f_{lm}(r)
30! {\bf Y}_{lm}^{l+1}(\hat{\bf r})\\
31! &+\sqrt{\frac{l}{2l+1}}\left(\frac{d}{dr}+\frac{l+1}{r}\right)f_{lm}(r)
32! {\bf Y}_{lm}^{l-1}(\hat{\bf r}),
33! \end{align*}
34! where the vector spherical harmonics are determined from Clebsch-Gordan
35! coefficients as follows:
36! $$ {\bf Y}_{lm}^{l'}(\hat{\bf r})=\sum_{m'\mu}
37! \begin{bmatrix} l' & 1 & l \\ m' & \mu & m \end{bmatrix}
38! Y_{lm}(\hat{\bf r})\hat{\bf e}^{\mu} $$
39! and the (contravariant) spherical unit vectors are given by
40! $$ \hat{\bf e}_{+1}=-\frac{\hat{\bf x}+i\hat{\bf y}}{\sqrt{2}},
41! \qquad\hat{\bf e}_0=\hat{\bf z},\qquad
42! \hat{\bf e}_{-1}=\frac{\hat{\bf x}-i\hat{\bf y}}{\sqrt{2}}. $$
43!
44! !REVISION HISTORY:
45! Rewritten May 2009 (JKD)
46! Modified, February 2020 (JKD)
47!EOP
48!BOC
49implicit none
50! arguments
51integer, intent(in) :: nr,nri
52real(8), intent(in) :: ri(nr),wcr(12,nr)
53complex(8), intent(in) :: zfmt(*)
54integer, intent(in) :: ld
55complex(8), intent(out) :: gzfmt(ld,3)
56! local variables
57integer nro,iro,ir,mu
58integer np,npi,i,i1,j
59integer l,m,lm,lm1
60! real constant 1/sqrt(2)
61real(8), parameter :: c1=0.7071067811865475244d0
62real(8) t1,t2,t3
63complex(8) z1
64! automatic arrays
65real(8) f1(nr),f2(nr),g1(nr),g2(nr)
66complex(8) drmt(ld)
67! external functions
68real(8), external :: clebgor
69nro=nr-nri
70iro=nri+1
71npi=lmmaxi*nri
72np=npi+lmmaxo*nro
73!----------------------------------------!
74! compute the radial derivatives !
75!----------------------------------------!
76do lm=1,lmmaxi
77 i=lm
78 do ir=1,nri
79 f1(ir)=dble(zfmt(i))
80 f2(ir)=aimag(zfmt(i))
81 i=i+lmmaxi
82 end do
83 do ir=iro,nr
84 f1(ir)=dble(zfmt(i))
85 f2(ir)=aimag(zfmt(i))
86 i=i+lmmaxo
87 end do
88 call splined(nr,wcr,f1,g1)
89 call splined(nr,wcr,f2,g2)
90 i=lm
91 do ir=1,nri
92 drmt(i)=cmplx(g1(ir),g2(ir),8)
93 i=i+lmmaxi
94 end do
95 do ir=iro,nr
96 drmt(i)=cmplx(g1(ir),g2(ir),8)
97 i=i+lmmaxo
98 end do
99end do
100do lm=lmmaxi+1,lmmaxo
101 i=npi+lm
102 do ir=iro,nr
103 f1(ir)=dble(zfmt(i))
104 f2(ir)=aimag(zfmt(i))
105 i=i+lmmaxo
106 end do
107 call splined(nro,wcr(1,iro),f1(iro),g1(iro))
108 call splined(nro,wcr(1,iro),f2(iro),g2(iro))
109 i=npi+lm
110 do ir=iro,nr
111 drmt(i)=cmplx(g1(ir),g2(ir),8)
112 i=i+lmmaxo
113 end do
114end do
115!-----------------------------------------------------!
116! compute the gradient in the spherical basis !
117!-----------------------------------------------------!
118! zero the gradient array
119gzfmt(1:np,1:3)=0.d0
120! inner part of muffin-tin
121lm=0
122do l=0,lmaxi
123 t1=-sqrt(dble(l+1)/dble(2*l+1))
124 if (l > 0) then
125 t2=sqrt(dble(l)/dble(2*l+1))
126 else
127 t2=0.d0
128 end if
129 do m=-l,l
130 lm=lm+1
131 j=1
132 do mu=-1,1
133 if (mu == 0) j=3
134 if (mu == 1) j=2
135 if (l+1 <= lmaxi) then
136! index to (l,m) is l*(l+1)+m+1, therefore index to (l+1,m-mu) is
137 lm1=(l+1)*(l+2)+(m-mu)+1
138 t3=t1*clebgor(l+1,1,l,m-mu,mu,m)
139 i=lm; i1=lm1
140 do ir=1,nri
141 gzfmt(i1,j)=gzfmt(i1,j)+t3*(drmt(i)-dble(l)*ri(ir)*zfmt(i))
142 i=i+lmmaxi; i1=i1+lmmaxi
143 end do
144 end if
145 if (abs(m-mu) <= l-1) then
146! index to (l-1,m-mu)
147 lm1=(l-1)*l+(m-mu)+1
148 t3=t2*clebgor(l-1,1,l,m-mu,mu,m)
149 i=lm; i1=lm1
150 do ir=1,nri
151 gzfmt(i1,j)=gzfmt(i1,j)+t3*(drmt(i)+dble(l+1)*ri(ir)*zfmt(i))
152 i=i+lmmaxi; i1=i1+lmmaxi
153 end do
154 end if
155 end do
156 end do
157end do
158! outer part of muffin-tin
159lm=0
160do l=0,lmaxo
161 t1=-sqrt(dble(l+1)/dble(2*l+1))
162 if (l > 0) then
163 t2=sqrt(dble(l)/dble(2*l+1))
164 else
165 t2=0.d0
166 end if
167 do m=-l,l
168 lm=lm+1
169 j=1
170 do mu=-1,1
171 if (mu == 0) j=3
172 if (mu == 1) j=2
173 if (l+1 <= lmaxo) then
174 lm1=(l+1)*(l+2)+(m-mu)+1
175 t3=t1*clebgor(l+1,1,l,m-mu,mu,m)
176 i=npi+lm; i1=npi+lm1
177 do ir=iro,nr
178 gzfmt(i1,j)=gzfmt(i1,j)+t3*(drmt(i)-dble(l)*ri(ir)*zfmt(i))
179 i=i+lmmaxo; i1=i1+lmmaxo
180 end do
181 end if
182 if (abs(m-mu) <= l-1) then
183 lm1=(l-1)*l+(m-mu)+1
184 t3=t2*clebgor(l-1,1,l,m-mu,mu,m)
185 i=npi+lm; i1=npi+lm1
186 do ir=iro,nr
187 gzfmt(i1,j)=gzfmt(i1,j)+t3*(drmt(i)+dble(l+1)*ri(ir)*zfmt(i))
188 i=i+lmmaxo; i1=i1+lmmaxo
189 end do
190 end if
191 end do
192 end do
193end do
194!---------------------------------------------------!
195! convert from spherical to Cartesian basis !
196!---------------------------------------------------!
197! note that the gradient transforms as a covariant vector, i.e. y -> -y
198i=0
199do ir=1,nri
200 do lm=1,lmmaxi
201 i=i+1
202 z1=gzfmt(i,1)
203 gzfmt(i,1)=c1*(z1-gzfmt(i,2))
204 gzfmt(i,2)=c1*(z1+gzfmt(i,2))*zmi
205 end do
206end do
207do ir=iro,nr
208 do lm=1,lmmaxo
209 i=i+1
210 z1=gzfmt(i,1)
211 gzfmt(i,1)=c1*(z1-gzfmt(i,2))
212 gzfmt(i,2)=c1*(z1+gzfmt(i,2))*zmi
213 end do
214end do
215end subroutine
216!EOC
217
real(8) function clebgor(j1, j2, j3, m1, m2, m3)
Definition clebgor.f90:10
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition gradzfmt.f90:10
complex(8), parameter zmi
Definition modmain.f90:1239
integer lmmaxi
Definition modmain.f90:207
integer lmaxo
Definition modmain.f90:201
integer lmaxi
Definition modmain.f90:205
integer lmmaxo
Definition modmain.f90:203
pure subroutine splined(n, wc, f, df)
Definition splined.f90:7