The Elk Code
genzvclmt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2010 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 genzvclmt(nrmt_,nrmti_,ld1,rl,wpr,ld2,zrhomt,zvclmt)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 integer, intent(in) :: nrmt_(nspecies),nrmti_(nspecies)
12 integer, intent(in) :: ld1
13 real(8), intent(in) :: rl(ld1,-lmaxo-1:lmaxo+2,nspecies)
14 real(8), intent(in) :: wpr(4,ld1,nspecies)
15 integer, intent(in) :: ld2
16 complex(8), intent(in) :: zrhomt(ld2,natmtot)
17 complex(8), intent(out) :: zvclmt(ld2,natmtot)
18 ! local variables
19 integer is,ias,nthd
20 call holdthd(natmtot,nthd)
21 !$OMP PARALLEL DO DEFAULT(SHARED) &
22 !$OMP PRIVATE(is) &
23 !$OMP SCHEDULE(DYNAMIC) &
24 !$OMP NUM_THREADS(nthd)
25 do ias=1,natmtot
26  is=idxis(ias)
27  call zpotclmt(nrmt_(is),nrmti_(is),ld1,rl(:,:,is),wpr(:,:,is),zrhomt(:,ias), &
28  zvclmt(:,ias))
29 end do
30 !$OMP END PARALLEL DO
31 call freethd(nthd)
32 end subroutine
33 
Definition: modomp.f90:6
subroutine genzvclmt(nrmt_, nrmti_, ld1, rl, wpr, ld2, zrhomt, zvclmt)
Definition: genzvclmt.f90:7
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 zpotclmt(nr, nri, ld, rl, wpr, zrhomt, zvclmt)
Definition: zpotclmt.f90:10