The Elk Code
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 
6 pure complex(8) function zcfmtinp(nr,nri,wr,cfmt1,cfmt2)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: nr,nri
11 real(8), intent(in) :: wr(nr)
12 complex(4), intent(in) :: cfmt1(*),cfmt2(*)
13 ! local variables
14 integer n,ir,i
15 ! compute the dot-products for each radial point and integrate over r
16 zcfmtinp=0.d0
17 if (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
27 else
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
34 end if
35 n=lmmaxo-1
36 do ir=nri+1,nr
37  zcfmtinp=zcfmtinp+wr(ir)*dot_product(cfmt1(i:i+n),cfmt2(i:i+n))
38  i=i+lmmaxo
39 end do
40 end function
41 
integer lmmaxo
Definition: modmain.f90:203
pure complex(8) function zcfmtinp(nr, nri, wr, cfmt1, cfmt2)
Definition: zcfmtinp.f90:7
integer lmmaxi
Definition: modmain.f90:207
integer lmaxi
Definition: modmain.f90:205