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