The Elk Code
wfmtsv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2021 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 wfmtsv(tsh,lrstp,is,ias,nst,idx,ngp,apwalm,evecfv,evecsv,ld,wfmt)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 logical, intent(in) :: tsh
12 integer, intent(in) :: lrstp,is,ias,nst,idx(*),ngp(nspnfv)
13 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv)
14 complex(8), intent(in) :: evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv)
15 integer, intent(in) :: ld
16 complex(8), intent(out) :: wfmt(ld,nspinor,nst)
17 ! local variables
18 logical tasv
19 integer io,ilo,ispn,jspn
20 integer nr,nri,nro,iro,np,npi
21 integer l,lm,lma,lmb,nm
22 integer n,p,i,j,k,nthd
23 complex(8) zq(2)
24 ! automatic arrays
25 complex(8) x(nstfv,nspnfv),y(nlmwf(is),nspinor,nst),z(2*lmaxo+1)
26 ! external functions
27 complex(8), external :: zdotu
28 iro=nrmti(is)+lrstp
29 if (lrstp == 1) then
30  nr=nrmt(is)
31  nri=nrmti(is)
32  np=npmt(is)
33  npi=npmti(is)
34 else
35  nr=nrcmt(is)
36  nri=nrcmti(is)
37  np=npcmt(is)
38  npi=npcmti(is)
39 end if
40 nro=nr-nri
41 ! de-phasing factor for spin-spirals
42 if (ssdph) then
43  zq(1)=zqss(ias)
44  zq(2)=conjg(zq(1))
45 end if
46 ! check if all the second-variational wavefunctions should be calculated
47 tasv=(idx(1) == 0)
48 call holdthd(nst,nthd)
49 !-----------------------!
50 ! APW functions !
51 !-----------------------!
52 !$OMP PARALLEL DEFAULT(SHARED) &
53 !$OMP PRIVATE(z,p,l,lma,lmb,io,lm) &
54 !$OMP PRIVATE(ispn,jspn,n,i,j,k,nm,ilo) &
55 !$OMP NUM_THREADS(nthd)
56 p=0
57 do l=0,lmaxo
58  lma=l**2+1; lmb=lma+2*l
59  do io=1,apword(l,is)
60  do lm=lma,lmb
61  p=p+1
62  if (tevecsv) then
63  do jspn=1,nspnfv
64  n=ngp(jspn)
65 !$OMP DO SCHEDULE(DYNAMIC)
66  do j=1,nstfv
67  x(j,jspn)=zdotu(n,evecfv(:,j,jspn),1,apwalm(:,io,lm,ias,jspn),1)
68  end do
69 !$OMP END DO
70  end do
71 ! loop only over required states
72 !$OMP DO SCHEDULE(DYNAMIC)
73  do j=1,nst
74 ! index to state in evecsv
75  k=merge(j,idx(j),tasv)
76  y(p,1,j)=zdotu(nstfv,evecsv(1,k),1,x,1)
77  if (spinpol) then
78  jspn=jspnfv(2)
79  y(p,2,j)=zdotu(nstfv,evecsv(nstfv+1,k),1,x(1,jspn),1)
80  end if
81  end do
82 !$OMP END DO
83  else
84 !$OMP DO SCHEDULE(DYNAMIC)
85  do j=1,nst
86  k=merge(j,idx(j),tasv)
87  y(p,1,j)=zdotu(ngp(1),evecfv(:,k,1),1,apwalm(:,io,lm,ias,1),1)
88  end do
89 !$OMP END DO
90  end if
91  end do
92  end do
93 end do
94 !$OMP DO SCHEDULE(DYNAMIC)
95 do j=1,nst
96  wfmt(1:np,1:nspinor,j)=0.d0
97  do ispn=1,nspinor
98  p=0
99  do l=0,lmaxo
100  nm=2*l+1
101  lma=l**2+1
102  do io=1,apword(l,is)
103  do lm=1,nm
104  p=p+1
105  z(lm)=y(p,ispn,j)
106  end do
107  if (ssdph) z(1:nm)=z(1:nm)*zq(ispn)
108  if (l <= lmaxi) then
109  call zfzrf(nm,nri,apwfr(1,1,io,l,ias),z,lmmaxi,wfmt(lma,ispn,j))
110  end if
111  call zfzrf(nm,nro,apwfr(iro,1,io,l,ias),z,lmmaxo,wfmt(npi+lma,ispn,j))
112  end do
113  end do
114  end do
115 end do
116 !$OMP END DO
117 !---------------------------------!
118 ! local-orbital functions !
119 !---------------------------------!
120 p=0
121 do ilo=1,nlorb(is)
122  l=lorbl(ilo,is)
123  do lm=l**2+1,(l+1)**2
124  p=p+1
125  i=idxlo(lm,ilo,ias)
126  if (tevecsv) then
127  do jspn=1,nspnfv
128  n=ngp(jspn)
129  x(1:nstfv,jspn)=evecfv(n+i,1:nstfv,jspn)
130  end do
131 !$OMP DO SCHEDULE(DYNAMIC)
132  do j=1,nst
133  k=merge(j,idx(j),tasv)
134  y(p,1,j)=zdotu(nstfv,evecsv(1,k),1,x,1)
135  if (spinpol) then
136  jspn=jspnfv(2)
137  y(p,2,j)=zdotu(nstfv,evecsv(nstfv+1,k),1,x(1,jspn),1)
138  end if
139  end do
140 !$OMP END DO
141  else
142  do j=1,nst
143  k=merge(j,idx(j),tasv)
144  y(p,1,j)=evecfv(ngp(1)+i,k,1)
145  end do
146  end if
147  end do
148 end do
149 !$OMP DO SCHEDULE(DYNAMIC)
150 do j=1,nst
151  do ispn=1,nspinor
152  p=0
153  do ilo=1,nlorb(is)
154  l=lorbl(ilo,is)
155  nm=2*l+1
156  lma=l**2+1
157  do lm=1,nm
158  p=p+1
159  z(lm)=y(p,ispn,j)
160  end do
161  if (ssdph) z(1:nm)=z(1:nm)*zq(ispn)
162  if (l <= lmaxi) then
163  call zfzrf(nm,nri,lofr(1,1,ilo,ias),z,lmmaxi,wfmt(lma,ispn,j))
164  end if
165  call zfzrf(nm,nro,lofr(iro,1,ilo,ias),z,lmmaxo,wfmt(npi+lma,ispn,j))
166  end do
167 ! convert to spherical coordinates if required
168  if (.not.tsh) call zbshtip(nr,nri,wfmt(:,ispn,j))
169  end do
170 end do
171 !$OMP END DO
172 !$OMP END PARALLEL
173 call freethd(nthd)
174 
175 contains
176 
177 pure subroutine zfzrf(m,n,rf,z,ld,zf)
178 implicit none
179 ! arguments
180 integer, intent(in) :: m,n
181 real(8), intent(in) :: rf(lrstp,n)
182 complex(8), intent(in) :: z(m)
183 integer, intent(in) :: ld
184 complex(8), intent(inout) :: zf(ld,n)
185 ! local variables
186 integer i
187 do i=1,m
188  zf(i,1:n)=zf(i,1:n)+z(i)*rf(1,1:n)
189 end do
190 end subroutine
191 
192 end subroutine
193 
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
logical spinpol
Definition: modmain.f90:228
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
logical tevecsv
Definition: modmain.f90:917
Definition: modomp.f90:6
complex(8), dimension(:), allocatable zqss
Definition: modmain.f90:287
subroutine zbshtip(nr, nri, zfmt)
Definition: zbshtip.f90:7
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
logical ssdph
Definition: modmain.f90:285
integer, dimension(maxspecies) npcmti
Definition: modmain.f90:214
integer lmmaxi
Definition: modmain.f90:207
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
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(maxspecies) npmti
Definition: modmain.f90:213
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:796
subroutine wfmtsv(tsh, lrstp, is, ias, nst, idx, ngp, apwalm, evecfv, evecsv, ld, wfmt)
Definition: wfmtsv.f90:7
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(2) jspnfv
Definition: modmain.f90:291
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150
integer lmaxi
Definition: modmain.f90:205