The Elk Code
zfmtinp.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: zfmtinp
8 ! !INTERFACE:
9 pure complex(8) function zfmtinp(nr,nri,wr,zfmt1,zfmt2)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! nr : number of radial mesh points (in,integer)
14 ! nri : number of points on the inner part of the muffin-tin (in,integer)
15 ! wr : weights for integration on radial mesh (in,real(nr))
16 ! zfmt1 : first complex muffin-tin function in spherical harmonics
17 ! (in,complex(*))
18 ! zfmt2 : second complex muffin-tin function (in,complex(*))
19 ! !DESCRIPTION:
20 ! Calculates the inner product of two complex fuctions in the muffin-tin. In
21 ! other words, given two complex functions of the form
22 ! $$ f({\bf r})=\sum_{l=0}^{l_{\rm max}}\sum_{m=-l}^{l}f_{lm}(r)Y_{lm}
23 ! (\hat{\bf r}), $$
24 ! the function returns
25 ! $$ I=\sum_{l=0}^{l_{\rm max}}\sum_{m=-l}^{l}\int f_{lm}^{1*}(r)
26 ! f_{lm}^2(r)r^2\,dr\;. $$
27 !
28 ! !REVISION HISTORY:
29 ! Created November 2003 (Sharma)
30 ! Modified, September 2013 (JKD)
31 ! Modified for packed functions, June 2016 (JKD)
32 !EOP
33 !BOC
34 implicit none
35 ! arguments
36 integer, intent(in) :: nr,nri
37 real(8), intent(in) :: wr(nr)
38 complex(8), intent(in) :: zfmt1(*),zfmt2(*)
39 ! local variables
40 integer n,ir,i
41 ! compute the dot-products for each radial point and integrate over r
42 zfmtinp=0.d0
43 if (lmaxi == 1) then
44  do ir=1,nri
45  i=4*(ir-1)+1
46  zfmtinp=zfmtinp+wr(ir) &
47  *(conjg(zfmt1(i))*zfmt2(i) &
48  +conjg(zfmt1(i+1))*zfmt2(i+1) &
49  +conjg(zfmt1(i+2))*zfmt2(i+2) &
50  +conjg(zfmt1(i+3))*zfmt2(i+3))
51  end do
52  i=4*nri+1
53 else
54  i=1
55  n=lmmaxi-1
56  do ir=1,nri
57  zfmtinp=zfmtinp+wr(ir)*dot_product(zfmt1(i:i+n),zfmt2(i:i+n))
58  i=i+lmmaxi
59  end do
60 end if
61 n=lmmaxo-1
62 do ir=nri+1,nr
63  zfmtinp=zfmtinp+wr(ir)*dot_product(zfmt1(i:i+n),zfmt2(i:i+n))
64  i=i+lmmaxo
65 end do
66 end function
67 !EOC
68 
integer lmmaxo
Definition: modmain.f90:203
pure complex(8) function zfmtinp(nr, nri, wr, zfmt1, zfmt2)
Definition: zfmtinp.f90:10
integer lmmaxi
Definition: modmain.f90:207
integer lmaxi
Definition: modmain.f90:205