The Elk Code
grad2rfmt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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: grad2rfmt
8 ! !INTERFACE:
9 subroutine grad2rfmt(nr,nri,ri,ri2,wcr,rfmt,g2rfmt)
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 ! ri : 1/r on the radial mesh (in,real(nr))
16 ! ri2 : 1/r^2 on the radial mesh (in,real(nr))
17 ! wcr : weights for spline coefficients on radial mesh (in,real(12,nr))
18 ! rfmt : real muffin-tin function (in,real(*))
19 ! g2rfmt : laplacian of the input function (out,real(*))
20 ! !DESCRIPTION:
21 ! Calculates the Laplacian of a real muffin-tin function. In other words,
22 ! given the real spherical harmonic expansion coefficients $f_{lm}(r)$ of a
23 ! function $f({\bf r})$, the routine returns
24 ! $$ F_{lm}(r)=\frac{1}{r}\frac{\partial^2}{\partial r^2}\big(rf_{lm}(r)\big)
25 ! -\frac{l(l+1)}{r^2}f_{lm}(r) $$
26 ! which yields
27 ! $$ \nabla^2 f({\bf r})=\sum_{lm} F_{lm}(r)R_{lm}(\hat{\bf r}), $$
28 ! where $R_{lm}$ is a real spherical harmonic function.
29 !
30 ! !REVISION HISTORY:
31 ! Created July 2009 (JKD)
32 !EOP
33 !BOC
34 implicit none
35 ! arguments
36 integer, intent(in) :: nr,nri
37 real(8), intent(in) :: ri(nr),ri2(nr),wcr(12,nr),rfmt(*)
38 real(8), intent(out) :: g2rfmt(*)
39 ! local variables
40 integer nro,iro,npi
41 integer l,lm,i1,j0,j1
42 real(8) t1
43 ! automatic arrays
44 real(8) fr(nr),cf(3,nr)
45 nro=nr-nri
46 iro=nri+1
47 npi=lmmaxi*nri
48 do l=0,lmaxi
49  t1=-dble(l*(l+1))
50  do lm=l**2+1,(l+1)**2
51 ! use a cubic spline to compute radial derivatives
52  i1=lmmaxi*(nri-1)+lm
53  j0=i1+lmmaxi
54  j1=lmmaxo*(nr-iro)+j0
55  fr(1:nri)=rfmt(lm:i1:lmmaxi)
56  fr(iro:nr)=rfmt(j0:j1:lmmaxo)
57  call splinew(nr,wcr,fr,cf)
58 ! apply Laplacian
59  g2rfmt(lm:i1:lmmaxi)=2.d0*(ri(1:nri)*cf(1,1:nri)+cf(2,1:nri)) &
60  +t1*ri2(1:nri)*rfmt(lm:i1:lmmaxi)
61  g2rfmt(j0:j1:lmmaxo)=2.d0*(ri(iro:nr)*cf(1,iro:nr)+cf(2,iro:nr)) &
62  +t1*ri2(iro:nr)*rfmt(j0:j1:lmmaxo)
63  end do
64 end do
65 do l=lmaxi+1,lmaxo
66  t1=-dble(l*(l+1))
67  do lm=l**2+1,(l+1)**2
68  j0=lmmaxi*nri+lm
69  j1=lmmaxo*(nr-iro)+j0
70  fr(iro:nr)=rfmt(j0:j1:lmmaxo)
71  call splinew(nro,wcr(:,iro),fr(iro),cf(1,iro))
72  g2rfmt(j0:j1:lmmaxo)=2.d0*(ri(iro:nr)*cf(1,iro:nr)+cf(2,iro:nr)) &
73  +t1*ri2(iro:nr)*rfmt(j0:j1:lmmaxo)
74  end do
75 end do
76 ! improve stability by smoothing the laplacian
77 call rfmtsm(msmgmt,nr,nri,g2rfmt)
78 end subroutine
79 !EOC
80 
integer lmmaxo
Definition: modmain.f90:203
pure subroutine splinew(n, wc, f, cf)
Definition: splinew.f90:7
integer msmgmt
Definition: modmain.f90:222
integer lmaxo
Definition: modmain.f90:201
subroutine rfmtsm(m, nr, nri, rfmt)
Definition: rfmtsm.f90:7
integer lmmaxi
Definition: modmain.f90:207
subroutine grad2rfmt(nr, nri, ri, ri2, wcr, rfmt, g2rfmt)
Definition: grad2rfmt.f90:10
integer lmaxi
Definition: modmain.f90:205