The Elk Code
 
Loading...
Searching...
No Matches
zftrf.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 Alexey I. Baranov.
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: zftrf
8! !INTERFACE:
9subroutine zftrf(npv,ivp,vpc,rfmt,rfir,zfp)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! npv : number of P-vectors (in,integer)
14! ivp : integer coordinates of the P-vectors (in,integer(3,npv))
15! vpc : P-vectors in Cartesian coordinates (in,real(3,npv))
16! rfmt : real muffin-tin function (in,real(npmtmax,natmtot))
17! rfir : real interstitial function (in,real(ngtot))
18! zfp : Fourier expansion coefficients of the real-space function
19! (out,complex(npv))
20! !DESCRIPTION:
21! Given a real function periodic in the unit cell, $f({\bf r})$, this routine
22! calculates its complex Fourier expansion coefficients:
23! $$ f({\bf P})=\frac{1}{\Omega}\int d^3r\,f({\bf r})\tilde{\Theta}({\bf r})
24! e^{-i{\bf P}\cdot{\bf r}}
25! +\frac{4\pi}{\Omega}\sum_{\alpha}e^{-i{\bf P}\cdot{\bf R}_{\alpha}}
26! \sum_{lm}(-i)^l Y_{lm}(\hat{\bf P})
27! \int_{0}^{R_{\alpha}}dr\,r^2 j_{l}(|{\bf P}|r)f_{lm}^{\alpha}(r), $$
28! where $\tilde{\Theta}$ is the smooth characteristic function of the
29! interstitial region, $\Omega$ is the unit cell volume and $R_{\alpha}$ is
30! the muffin-tin radius of atom $\alpha$.
31!
32! !REVISION HISTORY:
33! Created July 2010 (Alexey I. Baranov)
34! Modified, November 2010 (JKD)
35! Optimised, May 2024 (JKD)
36!EOP
37!BOC
38implicit none
39! arguments
40integer, intent(in) :: npv,ivp(3,npv)
41real(8), intent(in) :: vpc(3,npv),rfmt(npmtmax,natmtot),rfir(ngtot)
42complex(8), intent(out) :: zfp(npv)
43! local variables
44integer is,ia,ias,ip,ig,ifg
45integer nrc,nrci,irco,irc
46integer l,lm,n,i
47real(8) p,t0,t1
48complex(8) zsm,z1,z2
49! automatic arrays
50real(8) jl(0:lmaxo,nrcmtmax),rfmt1(npcmtmax)
51complex(8) ylm(lmmaxo)
52! allocatable arrays
53complex(8), allocatable :: zfft(:),zfmt(:,:)
54allocate(zfft(ngtot),zfmt(npcmtmax,natmtot))
55! zero the coefficients
56zfp(1:npv)=0.d0
57!---------------------------!
58! interstitial part !
59!---------------------------!
60! Fourier transform to G-space
61zfft(1:ngtot)=rfir(1:ngtot)
62call zfftifc(3,ngridg,-1,zfft)
63! find coefficients for all required input vectors
64do ip=1,npv
65 if ((ivp(1,ip) >= intgv(1,1)).and.(ivp(1,ip) <= intgv(2,1)).and. &
66 (ivp(2,ip) >= intgv(1,2)).and.(ivp(2,ip) <= intgv(2,2)).and. &
67 (ivp(3,ip) >= intgv(1,3)).and.(ivp(3,ip) <= intgv(2,3))) then
68 ig=ivgig(ivp(1,ip),ivp(2,ip),ivp(3,ip))
69 zfp(ip)=zfft(igfft(ig))
70 end if
71end do
72!-------------------------!
73! muffin-tin part !
74!-------------------------!
75! convert function from real to complex spherical harmonic expansion on coarse
76! radial mesh
77do ias=1,natmtot
78 is=idxis(ias)
79 nrc=nrcmt(is)
80 nrci=nrcmti(is)
81 call rfmtftoc(nrc,nrci,rfmt(:,ias),rfmt1)
82 call rtozfmt(nrc,nrci,rfmt1,zfmt(:,ias))
83end do
84! remove continuation of interstitial function into muffin-tin
85do ig=1,ngvec
86 ifg=igfft(ig)
87! conjugate spherical harmonics 4π iˡ Yₗₘ(G)*
88 call genylmv(.true.,lmaxo,vgc(:,ig),ylm)
89 ylm(1:lmmaxo)=conjg(ylm(1:lmmaxo))
90 do is=1,nspecies
91 nrc=nrcmt(is)
92 nrci=nrcmti(is)
93 irco=nrci+1
94! generate spherical Bessel functions
95 do irc=1,nrci
96 t1=gc(ig)*rcmt(irc,is)
97 call sbessel(lmaxi,t1,jl(:,irc))
98 end do
99 do irc=irco,nrc
100 t1=gc(ig)*rcmt(irc,is)
101 call sbessel(lmaxo,t1,jl(:,irc))
102 end do
103 do ia=1,natoms(is)
104 ias=idxas(ia,is)
105! structure factor
106 t1=vgc(1,ig)*atposc(1,ia,is) &
107 +vgc(2,ig)*atposc(2,ia,is) &
108 +vgc(3,ig)*atposc(3,ia,is)
109 z1=zfft(ifg)*cmplx(cos(t1),sin(t1),8)
110 i=1
111 do irc=1,nrci
112 do l=0,lmaxi
113 n=2*l
114 lm=l**2+1
115 z2=z1*jl(l,irc)
116 zfmt(i:i+n,ias)=zfmt(i:i+n,ias)-z2*ylm(lm:lm+n)
117 i=i+n+1
118 end do
119 end do
120 do irc=irco,nrc
121 do l=0,lmaxo
122 n=2*l
123 lm=l**2+1
124 z2=z1*jl(l,irc)
125 zfmt(i:i+n,ias)=zfmt(i:i+n,ias)-z2*ylm(lm:lm+n)
126 i=i+n+1
127 end do
128 end do
129 end do
130 end do
131end do
132t0=1.d0/omega
133! loop over input P-vectors
134do ip=1,npv
135! length of P-vector
136 p=sqrt(vpc(1,ip)**2+vpc(2,ip)**2+vpc(3,ip)**2)
137! spherical harmonics 4π (-i)ˡ Yₗₘ(P)
138 call genylmv(.true.,lmaxo,vpc(:,ip),ylm)
139 do is=1,nspecies
140 nrc=nrcmt(is)
141 nrci=nrcmti(is)
142 irco=nrci+1
143! generate spherical Bessel functions
144 do irc=1,nrci
145 t1=p*rcmt(irc,is)
146 call sbessel(lmaxi,t1,jl(:,irc))
147 end do
148 do irc=irco,nrc
149 t1=p*rcmt(irc,is)
150 call sbessel(lmaxo,t1,jl(:,irc))
151 end do
152 do ia=1,natoms(is)
153 ias=idxas(ia,is)
154 zsm=0.d0
155 i=1
156 do irc=1,nrci
157 z1=jl(0,irc)*zfmt(i,ias)*ylm(1)
158 i=i+1
159 do l=1,lmaxi
160 n=2*l
161 lm=l**2+1
162 z1=z1+jl(l,irc)*sum(zfmt(i:i+n,ias)*ylm(lm:lm+n))
163 i=i+n+1
164 end do
165 zsm=zsm+wr2cmt(irc,is)*z1
166 end do
167 do irc=irco,nrc
168 z1=jl(0,irc)*zfmt(i,ias)*ylm(1)
169 i=i+1
170 do l=1,lmaxo
171 n=2*l
172 lm=l**2+1
173 z1=z1+jl(l,irc)*sum(zfmt(i:i+n,ias)*ylm(lm:lm+n))
174 i=i+n+1
175 end do
176 zsm=zsm+wr2cmt(irc,is)*z1
177 end do
178! conjugate structure factor
179 t1=vpc(1,ip)*atposc(1,ia,is) &
180 +vpc(2,ip)*atposc(2,ia,is) &
181 +vpc(3,ip)*atposc(3,ia,is)
182 z1=t0*cmplx(cos(t1),-sin(t1),8)
183 zfp(ip)=zfp(ip)+z1*zsm
184 end do
185 end do
186end do
187deallocate(zfft,zfmt)
188end subroutine
189! EOC
190
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition genylmv.f90:10
integer, dimension(2, 3) intgv
Definition modmain.f90:394
integer, dimension(3) ngridg
Definition modmain.f90:386
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition modmain.f90:54
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(:,:), allocatable rcmt
Definition modmain.f90:177
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
real(8) omega
Definition modmain.f90:20
integer ngvec
Definition modmain.f90:396
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
integer lmaxi
Definition modmain.f90:205
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer nspecies
Definition modmain.f90:34
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
integer, dimension(:,:,:), allocatable ivgig
Definition modmain.f90:402
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition rfmtftoc.f90:7
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition rtozfmt.f90:7
subroutine sbessel(lmax, x, jl)
Definition sbessel.f90:10
subroutine zfftifc(nd, n, sgn, z)
subroutine zftrf(npv, ivp, vpc, rfmt, rfir, zfp)
Definition zftrf.f90:10