The Elk Code
 
Loading...
Searching...
No Matches
zcfinp.f90
Go to the documentation of this file.
1
2! Copyright (C) 2023 J. K. Dewhurst and S. Sharma.
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 zcfinp(cfmt1,cfir1,cfmt2,cfir2)
7use modmain
8use modomp
9implicit none
10! arguments
11complex(4), intent(in) :: cfmt1(npcmtmax,natmtot),cfir1(ngtc)
12complex(4), intent(in) :: cfmt2(npcmtmax,natmtot),cfir2(ngtc)
13! local variables
14integer is,ias,nthd
15complex(4) c1
16! external functions
17complex(4), external :: cdotc
18complex(8), external :: zcfmtinp
19zcfinp=0.d0
20call holdthd(natmtot+1,nthd)
21!$OMP PARALLEL DEFAULT(SHARED) &
22!$OMP PRIVATE(is,c1) REDUCTION(+:zcfinp) &
23!$OMP NUM_THREADS(nthd)
24!$OMP DO SCHEDULE(DYNAMIC)
25do ias=1,natmtot
26 is=idxis(ias)
27! muffin-tin contribution
28 zcfinp=zcfinp+zcfmtinp(nrcmt(is),nrcmti(is),wr2cmt(:,is),cfmt1(:,ias), &
29 cfmt2(:,ias))
30end do
31!$OMP END DO NOWAIT
32! interstitial contribution (requires that one of the functions has been
33! multiplied by the characteristic function)
34!$OMP SINGLE
35c1=cdotc(ngtc,cfir1,1,cfir2,1)
36zcfinp=zcfinp+c1*omega/dble(ngtc)
37!$OMP END SINGLE
38!$OMP END PARALLEL
39call freethd(nthd)
40end function
41
real(8) omega
Definition modmain.f90:20
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer ngtc
Definition modmain.f90:392
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
integer natmtot
Definition modmain.f90:40
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
pure complex(8) function zcfmtinp(nr, nri, wr, cfmt1, cfmt2)
Definition zcfmtinp.f90:7