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