The Elk Code
hmlfv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 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 hmlfv(nmatp,ngp,igpig,vgpc,apwalm,h)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 integer, intent(in) :: nmatp,ngp,igpig(ngkmax)
12 real(8), intent(in) :: vgpc(3,ngkmax)
13 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot)
14 complex(8), intent(out) :: h(nmatp,nmatp)
15 ! local variables
16 integer is,ias,j,nthd
17 ! interstitial contribution
18 call hmlistl(ngp,igpig,vgpc,nmatp,h)
19 ! APW-APW contribution
20 do ias=1,natmtot
21  is=idxis(ias)
22  call hmlaa(tefvr,is,ias,ngp,apwalm(:,:,:,ias),nmatp,h)
23 end do
24 ! zero the APW-lo and lo-lo part of the matrix
25 do j=ngp+1,nmatp
26  h(1:j,j)=0.d0
27 end do
28 call holdthd(natmtot,nthd)
29 !$OMP PARALLEL DO DEFAULT(SHARED) &
30 !$OMP PRIVATE(is) &
31 !$OMP SCHEDULE(DYNAMIC) &
32 !$OMP NUM_THREADS(nthd)
33 do ias=1,natmtot
34  is=idxis(ias)
35  call hmlalo(is,ias,ngp,apwalm(:,:,:,ias),nmatp,h)
36  call hmllolo(is,ias,ngp,nmatp,h)
37 end do
38 !$OMP END PARALLEL DO
39 call freethd(nthd)
40 end subroutine
41 
pure subroutine hmlistl(ngp, igpig, vgpc, ld, h)
Definition: hmlistl.f90:10
Definition: modomp.f90:6
subroutine hmlaa(thr, is, ias, ngp, apwalm, ld, h)
Definition: hmlaa.f90:10
logical tefvr
Definition: modmain.f90:870
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
pure subroutine hmllolo(is, ias, ngp, ld, h)
Definition: hmllolo.f90:7
subroutine hmlfv(nmatp, ngp, igpig, vgpc, apwalm, h)
Definition: hmlfv.f90:7
pure subroutine hmlalo(is, ias, ngp, apwalm, ld, h)
Definition: hmlalo.f90:7