The Elk Code
rfmtlm.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2016 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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: rfmtlm
8
! !INTERFACE:
9
pure subroutine
rfmtlm
(lm,nr,nri,rfmt,fr)
10
! !USES:
11
use
modmain
12
! lm : required (l,m) component (in,integer)
13
! nr : number of radial mesh points (in,integer)
14
! nri : number of points on inner part of muffin-tin (in,integer)
15
! rfmt : real muffin-tin function (in,real(npmtmax))
16
! fr : (l,m) component function on radial mesh (out,real(nr))
17
! !DESCRIPTION:
18
! Given a function expanded in real spherical harmonics
19
! $$ f({\bf r})=\sum_{l=0}^{l_{\rm max}}\sum_{m=-l}^l f_{lm}(r)
20
! R_{lm}(\hat{\bf r}), $$
21
! where $l_{\rm max}$ corresponds to {\tt lmaxi} and {\tt lmaxo} on the inner
22
! and outer parts of the muffin-tin, respectively, this routine returns a
23
! particular component function $f_{lm}$. See also {\tt genrlmv}.
24
!
25
! !REVISION HISTORY:
26
! Created April 2016 (JKD)
27
!EOP
28
!BOC
29
implicit none
30
! arguments
31
integer
,
intent(in)
:: lm,nr,nri
32
real(8)
,
intent(in)
:: rfmt(
npmtmax
)
33
real(8)
,
intent(out)
:: fr(nr)
34
! local variables
35
integer
iro,i0,i1
36
if
(lm >
lmmaxi
)
then
37
fr(1:nri)=0.d0
38
else
39
i1=
lmmaxi
*(nri-1)+lm
40
fr(1:nri)=rfmt(lm:i1:
lmmaxi
)
41
end if
42
iro=nri+1
43
if
(lm >
lmmaxo
)
then
44
fr(iro:nr)=0.d0
45
else
46
i0=
lmmaxi
*nri+lm
47
i1=
lmmaxo
*(nr-iro)+i0
48
fr(iro:nr)=rfmt(i0:i1:
lmmaxo
)
49
end if
50
end subroutine
51
!EOC
52
modmain::lmmaxo
integer lmmaxo
Definition:
modmain.f90:203
modmain
Definition:
modmain.f90:6
rfmtlm
pure subroutine rfmtlm(lm, nr, nri, rfmt, fr)
Definition:
rfmtlm.f90:10
modmain::lmmaxi
integer lmmaxi
Definition:
modmain.f90:207
modmain::npmtmax
integer npmtmax
Definition:
modmain.f90:216
rfmtlm.f90
Generated by
1.8.14