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