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