The Elk Code
dwfmtfv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2013 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 dwfmtfv(ias,ngp,ngpq,apwalmq,dapwalm,evecfv,devecfv,dwfmt)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 integer, intent(in) :: ias,ngp,ngpq
12 complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw)
13 complex(8), intent(in) :: dapwalm(ngkmax,apwordmax,lmmaxapw)
14 complex(8), intent(in) :: evecfv(nmatmax),devecfv(nmatmax)
15 complex(8), intent(out) :: dwfmt(*)
16 ! local variables
17 integer is,io,ilo
18 integer nrci,nrco,iro
19 integer l,lm,npci,i
20 complex(8) z
21 ! external functions
22 complex(8), external :: zdotu
23 is=idxis(ias)
24 iro=nrmti(is)+lradstp
25 nrci=nrcmti(is)
26 nrco=nrcmt(is)-nrci
27 npci=npcmti(is)
28 ! zero the wavefunction derivative
29 dwfmt(1:npcmt(is))=0.d0
30 !---------------------------------!
31 ! local-orbital functions !
32 !---------------------------------!
33 do ilo=1,nlorb(is)
34  l=lorbl(ilo,is)
35  do lm=l**2+1,(l+1)**2
36  i=npci+lm
37  z=devecfv(ngpq+idxlo(lm,ilo,ias))
38  if (l <= lmaxi) then
39  call zfzrf(nrci,lradstp,lofr(1,1,ilo,ias),lmmaxi,dwfmt(lm))
40  end if
41  call zfzrf(nrco,lradstp,lofr(iro,1,ilo,ias),lmmaxo,dwfmt(i))
42  end do
43 end do
44 !-----------------------!
45 ! APW functions !
46 !-----------------------!
47 do l=0,lmaxo
48  do lm=l**2+1,(l+1)**2
49  i=npci+lm
50  do io=1,apword(l,is)
51  z=zdotu(ngpq,devecfv,1,apwalmq(:,io,lm),1)
52  if (ias == iasph) then
53  z=z+zdotu(ngp,evecfv,1,dapwalm(:,io,lm),1)
54  end if
55  if (l <= lmaxi) then
56  call zfzrf(nrci,lradstp,apwfr(1,1,io,l,ias),lmmaxi,dwfmt(lm))
57  end if
58  call zfzrf(nrco,lradstp,apwfr(iro,1,io,l,ias),lmmaxo,dwfmt(i))
59  end do
60  end do
61 end do
62 
63 contains
64 
65 pure subroutine zfzrf(n,ld1,rf,ld2,zf)
66 implicit none
67 ! arguments
68 integer, intent(in) :: n
69 integer, intent(in) :: ld1
70 real(8), intent(in) :: rf(ld1,n)
71 integer, intent(in) :: ld2
72 complex(8), intent(inout) :: zf(ld2,n)
73 zf(1,1:n)=zf(1,1:n)+z*rf(1,1:n)
74 end subroutine
75 
76 end subroutine
77 
pure subroutine zfzrf(n, ld1, rf, ld2, zf)
Definition: dwfmtfv.f90:66
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:786
integer, dimension(:,:,:), allocatable idxlo
Definition: modmain.f90:855
real(8), dimension(:,:,:,:), allocatable lofr
Definition: modmain.f90:814
integer lmmaxo
Definition: modmain.f90:203
integer iasph
Definition: modphonon.f90:15
integer lmaxo
Definition: modmain.f90:201
subroutine dwfmtfv(ias, ngp, ngpq, apwalmq, dapwalm, evecfv, devecfv, dwfmt)
Definition: dwfmtfv.f90:7
integer lradstp
Definition: modmain.f90:171
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
integer, dimension(maxspecies) npcmti
Definition: modmain.f90:214
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
real(8), dimension(:,:,:,:,:), allocatable apwfr
Definition: modmain.f90:774
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:796
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer lmaxi
Definition: modmain.f90:205