The Elk Code
 
Loading...
Searching...
No Matches
zftcf.f90
Go to the documentation of this file.
1
2! Copyright (C) 2023 J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine zftcf(ngp,jlgpr,ylmgp,ld,sfacgp,cfmt,cfir,zfgp)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: ngp
11real(8), intent(in) :: jlgpr(njcmax,nspecies,ngp)
12complex(8), intent(in) :: ylmgp(lmmaxo,ngp)
13integer, intent(in) :: ld
14complex(8), intent(in) :: sfacgp(ld,natmtot)
15complex(4), intent(in) :: cfmt(npcmtmax,natmtot),cfir(ngtc)
16complex(8), intent(out) :: zfgp(ngp)
17! local variables
18integer is,ia,ias,ig
19integer nrc,nrci,irco,irc
20integer l,lm,n,i,j
21real(8) t0,y0
22complex(8) zsm,z1
23! automatic arrays
24complex(4) ylm(lmmaxo),cfft(ngtc)
25!-----------------------------------!
26! interstitial contribution !
27!-----------------------------------!
28! multiply by coarse characteristic function
29cfft(1:ngtc)=cfir(1:ngtc)*cfrc(1:ngtc)
30! Fourier transform to coarse G-grid
31call cfftifc(3,ngdgc,-1,cfft)
32zfgp(1:ngp)=cfft(igfc(1:ngp))
33!---------------------------------!
34! muffin-tin contribution !
35!---------------------------------!
36t0=1.d0/omega
37y0=t0*fourpi*y00
38do ig=1,ngp
39 ylm(1:lmmaxo)=t0*ylmgp(1:lmmaxo,ig)
40 do is=1,nspecies
41 nrc=nrcmt(is)
42 nrci=nrcmti(is)
43 irco=nrci+1
44 do ia=1,natoms(is)
45 ias=idxas(ia,is)
46 zsm=0.d0
47 i=1
48 j=1
49! inner part of muffin-tin, note that lmaxi >= 1
50 if (lmaxi == 1) then
51 do irc=1,nrci
52 zsm=zsm+wr2cmt(irc,is)*(jlgpr(j,is,ig)*y0*cfmt(i,ias) &
53 +jlgpr(j+1,is,ig)*(cfmt(i+1,ias)*ylm(2)+cfmt(i+2,ias)*ylm(3) &
54 +cfmt(i+3,ias)*ylm(4)))
55 i=i+4
56 j=j+2
57 end do
58 else
59 do irc=1,nrci
60 z1=jlgpr(j,is,ig)*y0*cfmt(i,ias)+jlgpr(j+1,is,ig) &
61 *(cfmt(i+1,ias)*ylm(2)+cfmt(i+2,ias)*ylm(3)+cfmt(i+3,ias)*ylm(4))
62 i=i+4
63 j=j+2
64 do l=2,lmaxi
65 n=2*l
66 lm=l**2+1
67 z1=z1+jlgpr(j,is,ig)*sum(cfmt(i:i+n,ias)*ylm(lm:lm+n))
68 i=i+n+1
69 j=j+1
70 end do
71 zsm=zsm+wr2cmt(irc,is)*z1
72 end do
73 end if
74! outer part of muffin-tin, note that lmaxo >= 3
75 do irc=irco,nrc
76 z1=jlgpr(j,is,ig)*y0*cfmt(i,ias)+jlgpr(j+1,is,ig) &
77 *(cfmt(i+1,ias)*ylm(2)+cfmt(i+2,ias)*ylm(3)+cfmt(i+3,ias)*ylm(4)) &
78 +jlgpr(j+2,is,ig) &
79 *(cfmt(i+4,ias)*ylm(5)+cfmt(i+5,ias)*ylm(6)+cfmt(i+6,ias)*ylm(7) &
80 +cfmt(i+7,ias)*ylm(8)+cfmt(i+8,ias)*ylm(9)) &
81 +jlgpr(j+3,is,ig) &
82 *(cfmt(i+9,ias)*ylm(10)+cfmt(i+10,ias)*ylm(11)+cfmt(i+11,ias)*ylm(12) &
83 +cfmt(i+12,ias)*ylm(13)+cfmt(i+13,ias)*ylm(14)+cfmt(i+14,ias)*ylm(15)&
84 +cfmt(i+15,ias)*ylm(16))
85 i=i+16
86 j=j+4
87 do l=4,lmaxo
88 n=2*l
89 lm=l**2+1
90 z1=z1+jlgpr(j,is,ig)*sum(cfmt(i:i+n,ias)*ylm(lm:lm+n))
91 i=i+n+1
92 j=j+1
93 end do
94 zsm=zsm+wr2cmt(irc,is)*z1
95 end do
96 zfgp(ig)=zfgp(ig)+conjg(sfacgp(ig,ias))*zsm
97 end do
98 end do
99end do
100end subroutine
101
subroutine cfftifc(nd, n, sgn, c)
real(8), parameter y00
Definition modmain.f90:1233
real(8), dimension(:), allocatable cfrc
Definition modmain.f90:438
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
real(8) omega
Definition modmain.f90:20
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
real(8), parameter fourpi
Definition modmain.f90:1231
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
integer lmaxo
Definition modmain.f90:201
integer lmaxi
Definition modmain.f90:205
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
subroutine zftcf(ngp, jlgpr, ylmgp, ld, sfacgp, cfmt, cfir, zfgp)
Definition zftcf.f90:7