The Elk Code
genexpmt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine genexpmt(ngp,jlgpr,ylmgp,ld,sfacgp,expmt)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: ngp
11 real(8), intent(in) :: jlgpr(njcmax,nspecies,ngp)
12 complex(8), intent(in) :: ylmgp(lmmaxo,ngp)
13 integer, intent(in) :: ld
14 complex(8), intent(in) :: sfacgp(ld,natmtot)
15 complex(8), intent(out) :: expmt(npcmtmax,natmtot,ngp)
16 ! local variables
17 integer ig,is,ia,ias
18 integer nrc,nrci,irc,npc
19 integer l,lma,lmb,i,j
20 real(8) t1
21 complex(8) z1
22 ! automatic arrays
23 complex(8) ylm(lmmaxo),zfmt(npcmtmax)
24 do ig=1,ngp
25  ylm(1:lmmaxo)=conjg(ylmgp(1:lmmaxo,ig))
26  do is=1,nspecies
27  nrc=nrcmt(is)
28  nrci=nrcmti(is)
29  npc=npcmt(is)
30  i=0
31  j=0
32  do irc=1,nrci
33  do l=0,lmaxi
34  j=j+1
35  t1=jlgpr(j,is,ig)
36  lma=l**2+1; lmb=lma+2*l
37  zfmt(i+lma:i+lmb)=t1*ylm(lma:lmb)
38  end do
39  i=i+lmmaxi
40  end do
41  do irc=nrci+1,nrc
42  do l=0,lmaxo
43  j=j+1
44  t1=jlgpr(j,is,ig)
45  lma=l**2+1; lmb=lma+2*l
46  zfmt(i+lma:i+lmb)=t1*ylm(lma:lmb)
47  end do
48  i=i+lmmaxo
49  end do
50 ! convert to spherical coordinates
51  call zbshtip(nrc,nrci,zfmt)
52 ! mutiply by phase factors and store for all atoms
53  do ia=1,natoms(is)
54  ias=idxas(ia,is)
55  z1=sfacgp(ig,ias)
56  expmt(1:npc,ias,ig)=z1*zfmt(1:npc)
57  end do
58  end do
59 end do
60 end subroutine
61 
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
integer lmaxo
Definition: modmain.f90:201
subroutine zbshtip(nr, nri, zfmt)
Definition: zbshtip.f90:7
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition: genexpmt.f90:7
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer lmmaxi
Definition: modmain.f90:207
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
integer lmaxi
Definition: modmain.f90:205