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