The Elk Code
Loading...
Searching...
No Matches
rfinpc.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2016 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
real(8)
function
rfinpc
(ld,rfmt1,rfir1,rfmt2,rfir2)
7
use
modmain
8
use
modomp
9
implicit none
10
! arguments
11
integer
,
intent(in)
:: ld
12
real
(8),
intent(in)
:: rfmt1(ld,
natmtot
),rfir1(
ngtot
)
13
real
(8),
intent(in)
:: rfmt2(ld,
natmtot
),rfir2(
ngtot
)
14
! local variables
15
integer
is,ias,nthd
16
! external functions
17
real
(8),
external
::
rfmtinp
18
! interstitial contribution
19
rfinpc
=sum(rfir1(1:
ngtot
)*rfir2(1:
ngtot
)*
cfunir
(1:
ngtot
))
20
rfinpc
=
rfinpc
*(
omega
/
ngtot
)
21
! muffin-tin contribution
22
call
holdthd
(
natmtot
,nthd)
23
!$OMP PARALLEL DO DEFAULT(SHARED) &
24
!$OMP PRIVATE(is) REDUCTION(+:rfinpc) &
25
!$OMP SCHEDULE(DYNAMIC) &
26
!$OMP NUM_THREADS(nthd)
27
do
ias=1,
natmtot
28
is=
idxis
(ias)
29
rfinpc
=
rfinpc
+
rfmtinp
(
nrcmt
(is),
nrcmti
(is),
wr2cmt
(:,is),rfmt1(:,ias), &
30
rfmt2(:,ias))
31
end do
32
!$OMP END PARALLEL DO
33
call
freethd
(nthd)
34
end function
35
modmain
Definition
modmain.f90:6
modmain::ngtot
integer ngtot
Definition
modmain.f90:390
modmain::omega
real(8) omega
Definition
modmain.f90:20
modmain::nrcmt
integer, dimension(maxspecies) nrcmt
Definition
modmain.f90:173
modmain::idxis
integer, dimension(maxatoms *maxspecies) idxis
Definition
modmain.f90:44
modmain::wr2cmt
real(8), dimension(:,:), allocatable wr2cmt
Definition
modmain.f90:189
modmain::natmtot
integer natmtot
Definition
modmain.f90:40
modmain::cfunir
real(8), dimension(:), allocatable cfunir
Definition
modmain.f90:436
modmain::nrcmti
integer, dimension(maxspecies) nrcmti
Definition
modmain.f90:211
modomp
Definition
modomp.f90:6
modomp::holdthd
subroutine holdthd(nloop, nthd)
Definition
modomp.f90:78
modomp::freethd
subroutine freethd(nthd)
Definition
modomp.f90:106
rfinpc
real(8) function rfinpc(ld, rfmt1, rfir1, rfmt2, rfir2)
Definition
rfinpc.f90:7
rfmtinp
pure real(8) function rfmtinp(nr, nri, wr, rfmt1, rfmt2)
Definition
rfmtinp.f90:10
rfinpc.f90
Generated by
1.9.8