The Elk Code
 
Loading...
Searching...
No Matches
rcfinp.f90
Go to the documentation of this file.
1
2! Copyright (C) 2014 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
6complex(8) function rcfinp(rfmt,rfir,cfmt,cfir)
7use modmain
8use modomp
9implicit none
10! arguments
11real(8), intent(in) :: rfmt(npcmtmax,natmtot),rfir(ngtot)
12complex(4), intent(in) :: cfmt(npcmtmax,natmtot),cfir(ngtot)
13! local variables
14integer is,ias,nthd
15! external functions
16complex(8), external :: rcfmtinp
17! interstitial contribution
18rcfinp=sum((cfunir(1:ngtot)*rfir(1:ngtot))*cfir(1:ngtot))
19rcfinp=rcfinp*(omega/ngtot)
20! muffin-tin contribution
21call holdthd(natmtot,nthd)
22!$OMP PARALLEL DO DEFAULT(SHARED) &
23!$OMP PRIVATE(is) REDUCTION(+:rcfinp) &
24!$OMP SCHEDULE(DYNAMIC) &
25!$OMP NUM_THREADS(nthd)
26do ias=1,natmtot
27 is=idxis(ias)
28 rcfinp=rcfinp+rcfmtinp(nrcmt(is),nrcmti(is),wr2cmt(:,is),rfmt(:,ias), &
29 cfmt(:,ias))
30end do
31!$OMP END PARALLEL DO
32call freethd(nthd)
33end function
34
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
pure complex(8) function rcfmtinp(nr, nri, wr, rfmt, cfmt)
Definition rcfmtinp.f90:7