The Elk Code
 
Loading...
Searching...
No Matches
genwfsv_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 genwfsv_sp(tsh,tgp,nst,idx,ngridg_,igfft_,ngp,igpig,apwalm,evecfv, &
7 evecsv,wfmt,ld,wfir)
8use modmain
9use modomp
10implicit none
11! arguments
12logical, intent(in) :: tsh,tgp
13integer, intent(in) :: nst,idx(*),ngridg_(3),igfft_(*)
14integer, intent(in) :: ngp(nspnfv),igpig(ngkmax,nspnfv)
15complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv)
16complex(8), intent(in) :: evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv)
17complex(4), intent(out) :: wfmt(npcmtmax,natmtot,nspinor,nst)
18integer, intent(in) :: ld
19complex(4), intent(out) :: wfir(ld,nspinor,nst)
20! local variables
21integer is,ias,ldmt,nthd
22if (nst < 1) return
23! muffin-tin wavefunction
24ldmt=npcmtmax*natmtot
25call holdthd(natmtot+1,nthd)
26!$OMP PARALLEL DEFAULT(SHARED) &
27!$OMP PRIVATE(ias,is) &
28!$OMP NUM_THREADS(nthd)
29!$OMP DO SCHEDULE(DYNAMIC)
30do ias=1,natmtot
31 is=idxis(ias)
32 call wfmtsv_sp(tsh,is,ias,nst,idx,ngp,apwalm,evecfv,evecsv,ldmt, &
33 wfmt(1,ias,1,1))
34end do
35!$OMP END DO NOWAIT
36! interstitial wavefunction
37!$OMP SINGLE
38call wfirsv_sp(tgp,nst,idx,ngridg_,igfft_,ngp,igpig,evecfv,evecsv,ld,wfir)
39!$OMP END SINGLE
40!$OMP END PARALLEL
41call freethd(nthd)
42end subroutine
43
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition genwfsv_sp.f90:8
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine wfirsv_sp(tgp, nst, idx, ngridg_, igfft_, ngp, igpig, evecfv, evecsv, ld, wfir)
Definition wfirsv_sp.f90:7
subroutine wfmtsv_sp(tsh, is, ias, nst, idx, ngp, apwalm, evecfv, evecsv, ld, wfmt)
Definition wfmtsv_sp.f90:7