The Elk Code
rfmtint.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2020 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 !BOP
7 ! !ROUTINE: rfmtint
8 ! !INTERFACE:
9 pure real(8) function rfmtint(nr,nri,wr,rfmt)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! nr : number of radial mesh points (in,integer)
14 ! nri : number of points on inner part of muffin-tin (in,integer)
15 ! wr : weights for integration on radial mesh (in,real(nr))
16 ! rfmt : real function inside muffin-tin (in,real(*))
17 ! !DESCRIPTION:
18 ! Computes the integral of a real muffin-tin function. In other words, if
19 ! $$ f({\bf r})=\sum_{l=0}^{l_{\rm max}}\sum_{m=-l}^{l}f_{lm}(r)R_{lm}
20 ! (\hat{\bf r}) $$
21 ! where $R_{lm}$ are the real spherical harmonics, then this routine returns
22 ! $$ I=4\pi Y_{00}\int f_{00}(r)r^2 dr\;. $$
23 !
24 ! !REVISION HISTORY:
25 ! Created July 2020 (JKD)
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 integer, intent(in) :: nr,nri
31 real(8), intent(in) :: wr(nr),rfmt(*)
32 ! local variables
33 integer iro,i0,i1
34 i1=lmmaxi*(nri-1)+1
35 rfmtint=sum(wr(1:nri)*rfmt(1:i1:lmmaxi))
36 iro=nri+1
37 i0=i1+lmmaxi
38 i1=lmmaxo*(nr-iro)+i0
39 rfmtint=rfmtint+sum(wr(iro:nr)*rfmt(i0:i1:lmmaxo))
41 end function
42 !EOC
43 
integer lmmaxo
Definition: modmain.f90:203
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition: rfmtint.f90:10
integer lmmaxi
Definition: modmain.f90:207
real(8), parameter fourpi
Definition: modmain.f90:1228
real(8), parameter y00
Definition: modmain.f90:1230