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