The Elk Code
 
Loading...
Searching...
No Matches
gendwfsv.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
6subroutine gendwfsv(tsh,tgp,nst,idx,ngp,ngpq,igpqig,apwalmq,dapwalm,evecfv, &
7 devecfv,evecsv,devecsv,dwfmt,ld,dwfir)
8use modmain
9implicit none
10! arguments
11logical, intent(in) :: tsh,tgp
12integer, intent(in) :: nst,idx(*)
13integer, intent(in) :: ngp(nspnfv),ngpq(nspnfv),igpqig(ngkmax,nspnfv)
14complex(8), intent(in) :: apwalmq(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv)
15complex(8), intent(in) :: dapwalm(ngkmax,apwordmax,lmmaxapw,nspnfv)
16complex(8), intent(in) :: evecfv(nmatmax,nstfv,nspnfv)
17complex(8), intent(in) :: devecfv(nmatmax,nstfv,nspnfv)
18complex(8), intent(in) :: evecsv(nstsv,nstsv),devecsv(nstsv,nstsv)
19complex(8), intent(out) :: dwfmt(npcmtmax,natmtot,nspinor,nst)
20integer, intent(in) :: ld
21complex(8), intent(out) :: dwfir(ld,nspinor,nst)
22! local variables
23integer ist,ispn,jspn
24integer is,ia,ias,nrc,nrci
25integer npc,igp,ifg,i,j,k
26real(8) t0
27complex(8) z1
28! automatic arrays
29logical done(nstfv),ddone(nstfv)
30! allocatable arrays
31complex(8), allocatable :: wfmt1(:,:),wfmt2(:),dwfmt1(:,:)
32!---------------------------------------------!
33! muffin-tin wavefunction derivatives !
34!---------------------------------------------!
35if (tevecsv) then
36 allocate(wfmt1(npcmtmax,nstfv),dwfmt1(npcmtmax,nstfv))
37end if
38if (.not.tsh) allocate(wfmt2(npcmtmax))
39do is=1,nspecies
40 nrc=nrcmt(is)
41 nrci=nrcmti(is)
42 npc=npcmt(is)
43 do ia=1,natoms(is)
44 ias=idxas(ia,is)
45 done(:)=.false.
46 do j=1,nst
47 k=idx(j)
48 if (tevecsv) then
49 i=0
50 do ispn=1,nspinor
51 jspn=jspnfv(ispn)
52 dwfmt(1:npc,ias,ispn,j)=0.d0
53 do ist=1,nstfv
54 i=i+1
55 z1=devecsv(i,k)
56!***** check if tq0 is needed here
57 if (abs(z1%re)+abs(z1%im) > epsocc) then
58 if (.not.done(ist)) then
59 if (tsh) then
60 call wfmtfv(ias,ngp(jspn),apwalmq(:,:,:,ias,jspn), &
61 evecfv(:,ist,jspn),wfmt1(:,ist))
62 else
63 call wfmtfv(ias,ngp(jspn),apwalmq(:,:,:,ias,jspn), &
64 evecfv(:,ist,jspn),wfmt2)
65 call zbsht(nrc,nrci,wfmt2,wfmt1(:,ist))
66 end if
67 done(ist)=.true.
68 end if
69 call zaxpy(npc,z1,wfmt1(:,ist),1,dwfmt(:,ias,ispn,j),1)
70 end if
71 z1=evecsv(i,k)
72 if (abs(z1%re)+abs(z1%im) > epsocc) then
73 if (.not.ddone(ist)) then
74 if (tsh) then
75 call dwfmtfv(ias,ngp(jspn),ngpq(jspn), &
76 apwalmq(:,:,:,ias,jspn),dapwalm(:,:,:,jspn), &
77 evecfv(:,ist,jspn),devecfv(:,ist,jspn),dwfmt1(:,ist))
78 else
79 call dwfmtfv(ias,ngp(jspn),ngpq(jspn), &
80 apwalmq(:,:,:,ias,jspn),dapwalm(:,:,:,jspn), &
81 evecfv(:,ist,jspn),devecfv(:,ist,jspn),wfmt2)
82 call zbsht(nrc,nrci,wfmt2,dwfmt1(:,ist))
83 end if
84 ddone(ist)=.true.
85 end if
86 call zaxpy(npc,z1,dwfmt1(:,ist),1,dwfmt(:,ias,ispn,j),1)
87 end if
88 end do
89 end do
90 else
91 if (tsh) then
92 call dwfmtfv(ias,ngp,ngpq,apwalmq(:,:,:,ias,1),dapwalm,evecfv(:,k,1),&
93 devecfv(:,k,1),dwfmt(:,ias,1,j))
94 else
95 call dwfmtfv(ias,ngp,ngpq,apwalmq(:,:,:,ias,1),dapwalm,evecfv(:,k,1),&
96 devecfv(:,k,1),wfmt2)
97 call zbsht(nrc,nrci,wfmt2,dwfmt(:,ias,1,j))
98 end if
99 end if
100 end do
101 end do
102end do
103if (tevecsv) deallocate(wfmt1,dwfmt1)
104if (.not.tsh) deallocate(wfmt2)
105!-----------------------------------------------!
106! interstitial wavefunction derivatives !
107!-----------------------------------------------!
108t0=1.d0/sqrt(omega)
109do j=1,nst
110 k=idx(j)
111 dwfir(1:ld,1:nspinor,j)=0.d0
112 if (tevecsv) then
113 i=0
114 do ispn=1,nspinor
115 jspn=jspnfv(ispn)
116 do ist=1,nstfv
117 i=i+1
118 z1=devecsv(i,k)
119!***** check if tq0 is needed here
120 if (abs(z1%re)+abs(z1%im) > epsocc) then
121 if (tgp) then
122 do igp=1,ngp(jspn)
123 dwfir(igp,ispn,j)=dwfir(igp,ispn,j)+z1*evecfv(igp,ist,jspn)
124 end do
125 else
126 z1=t0*z1
127 do igp=1,ngp(jspn)
128 ifg=igfc(igpqig(igp,jspn))
129 dwfir(ifg,ispn,j)=dwfir(ifg,ispn,j)+z1*evecfv(igp,ist,jspn)
130 end do
131 end if
132 end if
133 z1=evecsv(i,k)
134 if (abs(z1%re)+abs(z1%im) > epsocc) then
135 if (tgp) then
136 do igp=1,ngpq(jspn)
137 dwfir(igp,ispn,j)=dwfir(igp,ispn,j)+z1*devecfv(igp,ist,jspn)
138 end do
139 else
140 z1=t0*z1
141 do igp=1,ngpq(jspn)
142 ifg=igfc(igpqig(igp,jspn))
143 dwfir(ifg,ispn,j)=dwfir(ifg,ispn,j)+z1*devecfv(igp,ist,jspn)
144 end do
145 end if
146 end if
147 end do
148 end do
149 else
150 if (tgp) then
151 do igp=1,ngpq(1)
152 dwfir(igp,1,j)=devecfv(igp,k,1)
153 end do
154 else
155 do igp=1,ngpq(1)
156 ifg=igfc(igpqig(igp,1))
157 dwfir(ifg,1,j)=t0*devecfv(igp,k,1)
158 end do
159 end if
160 end if
161 if (.not.tgp) then
162 do ispn=1,nspinor
163 call zfftifc(3,ngdgc,1,dwfir(:,ispn,j))
164 end do
165 end if
166end do
167end subroutine
168
subroutine dwfmtfv(ias, ngp, ngpq, apwalmq, dapwalm, evecfv, devecfv, dwfmt)
Definition dwfmtfv.f90:7
subroutine gendwfsv(tsh, tgp, nst, idx, ngp, ngpq, igpqig, apwalmq, dapwalm, evecfv, devecfv, evecsv, devecsv, dwfmt, ld, dwfir)
Definition gendwfsv.f90:8
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
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer, dimension(2) jspnfv
Definition modmain.f90:291
real(8) epsocc
Definition modmain.f90:900
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer nspecies
Definition modmain.f90:34
logical tevecsv
Definition modmain.f90:920
subroutine wfmtfv(ias, ngp, apwalm, evecfv, wfmt)
Definition wfmtfv.f90:10
subroutine zbsht(nr, nri, zfmt1, zfmt2)
Definition zbsht.f90:10
subroutine zfftifc(nd, n, sgn, z)