The Elk Code
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:
9 subroutine gradzfmt(nr,nri,ri,wcr,zfmt,ld,gzfmt)
10 ! !USES:
11 use 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
49 implicit none
50 ! arguments
51 integer, intent(in) :: nr,nri
52 real(8), intent(in) :: ri(nr),wcr(12,nr)
53 complex(8), intent(in) :: zfmt(*)
54 integer, intent(in) :: ld
55 complex(8), intent(out) :: gzfmt(ld,3)
56 ! local variables
57 integer nro,iro,ir,mu
58 integer np,npi,i0,i1,i,j
59 integer l,m,lm,lm1
60 ! real constant 1/√2
61 real(8), parameter :: c1=0.7071067811865475244d0
62 real(8) t1,t2,t3
63 complex(8) z1
64 ! automatic arrays
65 complex(8) f(nr),df(nr),drmt(ld)
66 ! external functions
67 real(8), external :: clebgor
68 nro=nr-nri
69 iro=nri+1
70 npi=lmmaxi*nri
71 np=npi+lmmaxo*nro
72 !----------------------------------------!
73 ! compute the radial derivatives !
74 !----------------------------------------!
75 do lm=1,lmmaxi
76  i=lmmaxi*(nri-1)+lm
77  i0=i+lmmaxi
78  i1=lmmaxo*(nr-iro)+i0
79  f(1:nri)=zfmt(lm:i:lmmaxi)
80  f(iro:nr)=zfmt(i0:i1:lmmaxo)
81  call splined(nr,wcr,f,df)
82  drmt(lm:i:lmmaxi)=df(1:nri)
83  drmt(i0:i1:lmmaxo)=df(iro:nr)
84 end do
85 do lm=lmmaxi+1,lmmaxo
86  i0=lmmaxi*nri+lm
87  i1=lmmaxo*(nr-iro)+i0
88  f(iro:nr)=zfmt(i0:i1:lmmaxo)
89  call splined(nro,wcr(1,iro),f(iro),df(iro))
90  drmt(i0:i1:lmmaxo)=df(iro:nr)
91 end do
92 !-----------------------------------------------------!
93 ! compute the gradient in the spherical basis !
94 !-----------------------------------------------------!
95 ! zero the gradient array
96 gzfmt(1:np,1:3)=0.d0
97 ! inner part of muffin-tin
98 lm=0
99 do l=0,lmaxi
100  t1=-sqrt(dble(l+1)/dble(2*l+1))
101  t2=merge(sqrt(dble(l)/dble(2*l+1)),0.d0,l > 0)
102  do m=-l,l
103  lm=lm+1
104  j=1
105  do mu=-1,1
106  if (mu == 0) j=3
107  if (mu == 1) j=2
108  if (l+1 <= lmaxi) then
109 ! index to (l,m) is l*(l+1)+m+1, therefore index to (l+1,m-mu) is
110  lm1=(l+1)*(l+2)+(m-mu)+1
111  t3=t1*clebgor(l+1,1,l,m-mu,mu,m)
112  i=lm; i1=lm1
113  do ir=1,nri
114  gzfmt(i1,j)=gzfmt(i1,j)+t3*(drmt(i)-dble(l)*ri(ir)*zfmt(i))
115  i=i+lmmaxi; i1=i1+lmmaxi
116  end do
117  end if
118  if (abs(m-mu) <= l-1) then
119 ! index to (l-1,m-mu)
120  lm1=(l-1)*l+(m-mu)+1
121  t3=t2*clebgor(l-1,1,l,m-mu,mu,m)
122  i=lm; i1=lm1
123  do ir=1,nri
124  gzfmt(i1,j)=gzfmt(i1,j)+t3*(drmt(i)+dble(l+1)*ri(ir)*zfmt(i))
125  i=i+lmmaxi; i1=i1+lmmaxi
126  end do
127  end if
128  end do
129  end do
130 end do
131 ! outer part of muffin-tin
132 lm=0
133 do l=0,lmaxo
134  t1=-sqrt(dble(l+1)/dble(2*l+1))
135  t2=merge(sqrt(dble(l)/dble(2*l+1)),0.d0,l > 0)
136  do m=-l,l
137  lm=lm+1
138  j=1
139  do mu=-1,1
140  if (mu == 0) j=3
141  if (mu == 1) j=2
142  if (l+1 <= lmaxo) then
143  lm1=(l+1)*(l+2)+(m-mu)+1
144  t3=t1*clebgor(l+1,1,l,m-mu,mu,m)
145  i=npi+lm; i1=npi+lm1
146  do ir=iro,nr
147  gzfmt(i1,j)=gzfmt(i1,j)+t3*(drmt(i)-dble(l)*ri(ir)*zfmt(i))
148  i=i+lmmaxo; i1=i1+lmmaxo
149  end do
150  end if
151  if (abs(m-mu) <= l-1) then
152  lm1=(l-1)*l+(m-mu)+1
153  t3=t2*clebgor(l-1,1,l,m-mu,mu,m)
154  i=npi+lm; i1=npi+lm1
155  do ir=iro,nr
156  gzfmt(i1,j)=gzfmt(i1,j)+t3*(drmt(i)+dble(l+1)*ri(ir)*zfmt(i))
157  i=i+lmmaxo; i1=i1+lmmaxo
158  end do
159  end if
160  end do
161  end do
162 end do
163 !---------------------------------------------------!
164 ! convert from spherical to Cartesian basis !
165 !---------------------------------------------------!
166 ! note that the gradient transforms as a covariant vector, i.e. y -> -y
167 i=0
168 do ir=1,nri
169  do lm=1,lmmaxi
170  i=i+1
171  z1=gzfmt(i,1)
172  gzfmt(i,1)=c1*(z1-gzfmt(i,2))
173  z1=c1*(z1+gzfmt(i,2))
174  gzfmt(i,2)=cmplx(z1%im,-z1%re,8)
175  end do
176 end do
177 do ir=iro,nr
178  do lm=1,lmmaxo
179  i=i+1
180  z1=gzfmt(i,1)
181  gzfmt(i,1)=c1*(z1-gzfmt(i,2))
182  z1=c1*(z1+gzfmt(i,2))
183  gzfmt(i,2)=cmplx(z1%im,-z1%re,8)
184  end do
185 end do
186 
187 contains
188 
189 pure subroutine splined(n,wc,f,df)
190 implicit none
191 ! arguments
192 integer, intent(in) :: n
193 real(8), intent(in) :: wc(12,n)
194 complex(8), intent(in) :: f(n)
195 complex(8), intent(out) :: df(n)
196 ! local variables
197 integer i
198 df(1)=wc(1,1)*f(1)+wc(2,1)*f(2)+wc(3,1)*f(3)+wc(4,1)*f(4)
199 df(2)=wc(1,2)*f(1)+wc(2,2)*f(2)+wc(3,2)*f(3)+wc(4,2)*f(4)
200 do i=3,n-2
201  df(i)=wc(1,i)*f(i-1)+wc(2,i)*f(i)+wc(3,i)*f(i+1)+wc(4,i)*f(i+2)
202 end do
203 i=n-1
204 df(i)=wc(1,i)*f(n-3)+wc(2,i)*f(n-2)+wc(3,i)*f(n-1)+wc(4,i)*f(n)
205 df(n)=wc(1,n)*f(n-3)+wc(2,n)*f(n-2)+wc(3,n)*f(n-1)+wc(4,n)*f(n)
206 end subroutine
207 
208 end subroutine
209 !EOC
210 
integer lmmaxo
Definition: modmain.f90:203
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition: gradzfmt.f90:10
integer lmaxo
Definition: modmain.f90:201
integer lmmaxi
Definition: modmain.f90:207
pure subroutine splined(n, wc, f, df)
Definition: gensocfr.f90:66
integer lmaxi
Definition: modmain.f90:205