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