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