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