File: zftrf.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (196 lines) | stat: -rw-r--r-- 5,523 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196

! Copyright (C) 2010 Alexey I. Baranov.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: zftrf
! !INTERFACE:
subroutine zftrf(npv,ivp,vpc,rfmt,rfir,zfp)
! !USES:
use modmain
! !INPUT/OUTPUT PARAMETERS:
!   npv  : number of P-vectors (in,integer)
!   ivp  : integer coordinates of the P-vectors (in,integer(3,npv))
!   vpc  : P-vectors in Cartesian coordinates (in,real(3,npv))
!   rfmt : real muffin-tin function (in,real(npmtmax,natmtot))
!   rfir : real interstitial function (in,real(ngtot))
!   zfp  : Fourier expansion coefficients of the real-space function
!          (out,complex(npv))
! !DESCRIPTION:
!   Given a real function periodic in the unit cell, $f({\bf r})$, this routine
!   calculates its complex Fourier expansion coefficients:
!   $$ f({\bf P})=\frac{1}{\Omega}\int d^3r\,f({\bf r})\tilde{\Theta}({\bf r})
!    e^{-i{\bf P}\cdot{\bf r}}
!    +\frac{4\pi}{\Omega}\sum_{\alpha}e^{-i{\bf P}\cdot{\bf R}_{\alpha}}
!    \sum_{lm}(-i)^l Y_{lm}(\hat{\bf P})
!    \int_{0}^{R_{\alpha}}dr\,r^2 j_{l}(|{\bf P}|r)f_{lm}^{\alpha}(r), $$
!   where $\tilde{\Theta}$ is the smooth characteristic function of the
!   interstitial region, $\Omega$ is the unit cell volume and $R_{\alpha}$ is
!   the muffin-tin radius of atom $\alpha$.
!
! !REVISION HISTORY:
!   Created July 2010 (Alexey I. Baranov)
!   Modified, November 2010 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: npv,ivp(3,npv)
real(8), intent(in) :: vpc(3,npv)
real(8), intent(in) :: rfmt(npmtmax,natmtot),rfir(ngtot)
complex(8), intent(out) :: zfp(npv)
! local variables
integer is,ia,ias
integer nrc,nrci,irco,irc
integer lmax,l,m,lm,npci,i
integer ip,ig,ifg
real(8) p,tp(2),t0,t1,t2,t3
complex(8) zsum1,zsum2
complex(8) z1,z2,z3
! automatic arrays
real(8) jl(0:lmaxo,nrcmtmax)
real(8) fr1(nrcmtmax),fr2(nrcmtmax),rfmt1(npcmtmax)
complex(8) ylm(lmmaxo)
! allocatable arrays
complex(8), allocatable :: zfft(:),zfmt(:,:)
! external functions
real(8) splint
external splint
allocate(zfft(ngtot),zfmt(npcmtmax,natmtot))
! zero the coefficients
zfp(:)=0.d0
!---------------------------!
!     interstitial part     !
!---------------------------!
! Fourier transform to G-space
zfft(:)=rfir(:)
call zfftifc(3,ngridg,-1,zfft)
! find coefficients for all required input vectors
do ip=1,npv
  if ((ivp(1,ip).ge.intgv(1,1)).and.(ivp(1,ip).le.intgv(2,1)).and. &
      (ivp(2,ip).ge.intgv(1,2)).and.(ivp(2,ip).le.intgv(2,2)).and. &
      (ivp(3,ip).ge.intgv(1,3)).and.(ivp(3,ip).le.intgv(2,3))) then
    ig=ivgig(ivp(1,ip),ivp(2,ip),ivp(3,ip))
    zfp(ip)=zfft(igfft(ig))
  end if
end do
!-------------------------!
!     muffin-tin part     !
!-------------------------!
! convert function from real to complex spherical harmonic expansion on coarse
! radial mesh
do ias=1,natmtot
  is=idxis(ias)
  call rfmtftoc(nrmt(is),nrmti(is),rfmt(:,ias),rfmt1)
  call rtozfmt(nrcmt(is),nrcmti(is),rfmt1,zfmt(:,ias))
end do
! remove continuation of interstitial function into muffin-tin
do ig=1,ngtot
  ifg=igfft(ig)
! spherical coordinates of G
  call sphcrd(vgc(:,ig),t1,tp)
! conjugate spherical harmonics Y_lm*(G)
  call genylm(lmaxo,tp,ylm)
  ylm(:)=conjg(ylm(:))
  do is=1,nspecies
    nrc=nrcmt(is)
    nrci=nrcmti(is)
    irco=nrci+1
    npci=npcmti(is)
! generate spherical Bessel functions
    lmax=lmaxi
    do irc=1,nrc
      t1=gc(ig)*rcmt(irc,is)
      call sbessel(lmax,t1,jl(:,irc))
      if (irc.eq.nrci) lmax=lmaxo
    end do
    do ia=1,natoms(is)
      ias=idxas(ia,is)
! structure factor
      t1=dot_product(vgc(:,ig),atposc(:,ia,is))
      z1=fourpi*zfft(ifg)*cmplx(cos(t1),sin(t1),8)
      lm=0
      do l=0,lmaxi
        z2=z1*zil(l)
        do m=-l,l
          lm=lm+1
          z3=z2*ylm(lm)
          i=lm
          do irc=1,nrci
            zfmt(i,ias)=zfmt(i,ias)-z3*jl(l,irc)
            i=i+lmmaxi
          end do
        end do
      end do
      lm=0
      do l=0,lmaxo
        z2=z1*zil(l)
        do m=-l,l
          lm=lm+1
          z3=z2*ylm(lm)
          i=npci+lm
          do irc=irco,nrc
            zfmt(i,ias)=zfmt(i,ias)-z3*jl(l,irc)
            i=i+lmmaxo
          end do
        end do
      end do
    end do
  end do
end do
t0=fourpi/omega
! loop over input P-vectors
do ip=1,npv
! spherical coordinates of P
  call sphcrd(vpc(:,ip),p,tp)
! spherical harmonics Y_lm(P)
  call genylm(lmaxo,tp,ylm)
  do is=1,nspecies
    nrc=nrcmt(is)
    nrci=nrcmti(is)
! generate spherical Bessel functions
    lmax=lmaxi
    do irc=1,nrc
      t1=p*rcmt(irc,is)
      call sbessel(lmax,t1,jl(:,irc))
      if (irc.eq.nrci) lmax=lmaxo
    end do
    do ia=1,natoms(is)
      ias=idxas(ia,is)
      lmax=lmaxi
      i=0
      do irc=1,nrc
        i=i+1
        zsum1=jl(0,irc)*zfmt(i,ias)*ylm(1)
        lm=1
        do l=1,lmax
          lm=lm+1
          i=i+1
          zsum2=zfmt(i,ias)*ylm(lm)
          do m=1-l,l
            lm=lm+1
            i=i+1
            zsum2=zsum2+zfmt(i,ias)*ylm(lm)
          end do
          zsum1=zsum1+jl(l,irc)*zilc(l)*zsum2
        end do
        zsum1=zsum1*r2cmt(irc,is)
        fr1(irc)=dble(zsum1)
        fr2(irc)=aimag(zsum1)
        if (irc.eq.nrci) lmax=lmaxo
      end do
      t1=splint(nrc,rcmt(:,is),fr1)
      t2=splint(nrc,rcmt(:,is),fr2)
! conjugate structure factor
      t3=-dot_product(vpc(:,ip),atposc(:,ia,is))
      z1=t0*cmplx(cos(t3),sin(t3),8)
      zfp(ip)=zfp(ip)+z1*cmplx(t1,t2,8)
    end do
  end do
end do
deallocate(zfft,zfmt)
return
end subroutine
! EOC