The Elk Code
wfmtfv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: wfmtfv
8 ! !INTERFACE:
9 subroutine wfmtfv(ias,ngp,apwalm,evecfv,wfmt)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! ias : joint atom and species number (in,integer)
14 ! ngp : number of G+p-vectors (in,integer)
15 ! apwalm : APW matching coefficients (in,complex(ngkmax,apwordmax,lmmaxapw))
16 ! evecfv : first-variational eigenvector (in,complex(nmatmax))
17 ! wfmt : complex muffin-tin wavefunction passed in as real array
18 ! (out,real(2,*))
19 ! !DESCRIPTION:
20 ! Calculates the first-variational wavefunction in the muffin-tin in terms of
21 ! a spherical harmonic expansion. For atom $\alpha$ and a particular $k$-point
22 ! ${\bf p}$, the $r$-dependent $(l,m)$-coefficients of the wavefunction for
23 ! the $i$th state are given by
24 ! $$ \Phi^{i{\bf p}}_{\alpha lm}(r)=\sum_{\bf G}b^{i{\bf p}}_{\bf G}
25 ! \sum_{j=1}^{M^{\alpha}_l}A^{\alpha}_{jlm}({\bf G+p})u^{\alpha}_{jl}(r)
26 ! +\sum_{j=1}^{N^{\alpha}}b^{i{\bf p}}_{(\alpha,j,m)}v^{\alpha}_j(r)
27 ! \delta_{l,l_j}, $$
28 ! where $b^{i{\bf p}}$ is the $i$th eigenvector returned from routine
29 ! {\tt eveqn}; $A^{\alpha}_{jlm}({\bf G+p})$ is the matching coefficient;
30 ! $M^{\alpha}_l$ is the order of the APW; $u^{\alpha}_{jl}$ is the APW radial
31 ! function; $N^{\alpha}$ is the number of local-orbitals; $v^{\alpha}_j$ is
32 ! the $j$th local-orbital radial function; and $(\alpha,j,m)$ is a compound
33 ! index for the location of the local-orbital in the eigenvector. See routines
34 ! {\tt genapwfr}, {\tt genlofr}, {\tt match} and {\tt eveqn}.
35 !
36 ! !REVISION HISTORY:
37 ! Created April 2003 (JKD)
38 ! Fixed description, October 2004 (C. Brouder)
39 ! Removed argument ist, November 2006 (JKD)
40 ! Changed arguments and optimised, December 2014 (JKD)
41 !EOP
42 !BOC
43 implicit none
44 ! arguments
45 integer, intent(in) :: ias,ngp
46 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw),evecfv(nmatmax)
47 complex(8), intent(out) :: wfmt(*)
48 ! local variables
49 integer is,io,ilo,ld
50 integer nrci,nrco,iro
51 integer l,lma,nm,npci
52 ! automatic arrays
53 complex(8) z(2*lmaxo+1)
54 is=idxis(ias)
55 iro=nrmti(is)+lradstp
56 nrci=nrcmti(is)
57 nrco=nrcmt(is)-nrci
58 npci=npcmti(is)
59 ! zero the wavefunction
60 wfmt(1:npcmt(is))=0.d0
61 !---------------------------------!
62 ! local-orbital functions !
63 !---------------------------------!
64 do ilo=1,nlorb(is)
65  l=lorbl(ilo,is)
66  nm=2*l+1
67  lma=l**2+1
68  z(1:nm)=evecfv(ngp+idxlo(lma:lma+2*l,ilo,ias))
69  if (l <= lmaxi) call zfzrf(nm,nrci,lofr(1,1,ilo,ias),z,lmmaxi,wfmt(lma))
70  call zfzrf(nm,nrco,lofr(iro,1,ilo,ias),z,lmmaxo,wfmt(npci+lma))
71 end do
72 !-----------------------!
73 ! APW functions !
74 !-----------------------!
75 ld=ngkmax*apwordmax
76 do l=0,lmaxo
77  nm=2*l+1
78  lma=l**2+1
79  do io=1,apword(l,is)
80  call zgemv('T',ngp,nm,zone,apwalm(:,io,lma),ld,evecfv,1,zzero,z,1)
81  if (l <= lmaxi) call zfzrf(nm,nrci,apwfr(1,1,io,l,ias),z,lmmaxi,wfmt(lma))
82  call zfzrf(nm,nrco,apwfr(iro,1,io,l,ias),z,lmmaxo,wfmt(npci+lma))
83  end do
84 end do
85 
86 contains
87 
88 pure subroutine zfzrf(m,n,rf,z,ld,zf)
89 implicit none
90 ! arguments
91 integer, intent(in) :: m,n
92 real(8), intent(in) :: rf(lradstp,n)
93 complex(8), intent(in) :: z(m)
94 integer, intent(in) :: ld
95 complex(8), intent(inout) :: zf(ld,n)
96 ! local variables
97 integer i
98 do i=1,m
99  zf(i,1:n)=zf(i,1:n)+z(i)*rf(1,1:n)
100 end do
101 end subroutine
102 
103 end subroutine
104 !EOC
105 
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:853
real(8), dimension(:,:,:,:), allocatable lofr
Definition: modmain.f90:814
integer lmmaxo
Definition: modmain.f90:203
complex(8), parameter zone
Definition: modmain.f90:1234
integer lradstp
Definition: modmain.f90:171
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
complex(8), parameter zzero
Definition: modmain.f90:1234
subroutine wfmtfv(ias, ngp, apwalm, evecfv, wfmt)
Definition: wfmtfv.f90:10
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