The Elk Code
 
Loading...
Searching...
No Matches
zcfmtinp.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
6pure complex(8) function zcfmtinp(nr,nri,wr,cfmt1,cfmt2)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: nr,nri
11real(8), intent(in) :: wr(nr)
12complex(4), intent(in) :: cfmt1(*),cfmt2(*)
13! local variables
14integer n,ir,i
15! compute the dot-products for each radial point and integrate over r
16zcfmtinp=0.d0
17if (lmaxi == 1) then
18 do ir=1,nri
19 i=4*(ir-1)+1
20 zcfmtinp=zcfmtinp+wr(ir) &
21 *(conjg(cfmt1(i))*cfmt2(i) &
22 +conjg(cfmt1(i+1))*cfmt2(i+1) &
23 +conjg(cfmt1(i+2))*cfmt2(i+2) &
24 +conjg(cfmt1(i+3))*cfmt2(i+3))
25 end do
26 i=4*nri+1
27else
28 i=1
29 n=lmmaxi-1
30 do ir=1,nri
31 zcfmtinp=zcfmtinp+wr(ir)*dot_product(cfmt1(i:i+n),cfmt2(i:i+n))
32 i=i+lmmaxi
33 end do
34end if
35n=lmmaxo-1
36do ir=nri+1,nr
37 zcfmtinp=zcfmtinp+wr(ir)*dot_product(cfmt1(i:i+n),cfmt2(i:i+n))
38 i=i+lmmaxo
39end do
40end function
41
integer lmmaxi
Definition modmain.f90:207
integer lmaxi
Definition modmain.f90:205
integer lmmaxo
Definition modmain.f90:203
pure complex(8) function zcfmtinp(nr, nri, wr, cfmt1, cfmt2)
Definition zcfmtinp.f90:7