The Elk Code
gradrfmt.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 Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: gradrfmt
8 ! !INTERFACE:
9 subroutine gradrfmt(nr,nri,ri,wcr,rfmt,ld,grfmt)
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 ! ri : 1/r on the radial mesh (in,real(nr))
16 ! wcr : weights for spline coefficients on radial mesh (in,real(12,nr))
17 ! rfmt : real muffin-tin function (in,real(*))
18 ! ld : leading dimension (in,integer)
19 ! grfmt : gradient of rfmt (out,real(ld,3))
20 ! !DESCRIPTION:
21 ! Calculates the gradient of a real muffin-tin function. In other words, given
22 ! the real spherical harmonic expansion coefficients $f_{lm}(r)$ of a function
23 ! $f({\bf r})$, the routine returns ${\bf F}_{lm}$ where
24 ! $$ \sum_{lm}{\bf F}_{lm}(r)R_{lm}(\hat{\bf r})=\nabla f({\bf r}), $$
25 ! and $R_{lm}$ is a real spherical harmonic function. This is done by first
26 ! converting the function to a complex spherical harmonic expansion and then
27 ! using the routine {\tt gradzfmt}. See routine {\tt genrlm}.
28 !
29 ! !REVISION HISTORY:
30 ! Created August 2003 (JKD)
31 !EOP
32 !BOC
33 implicit none
34 ! arguments
35 integer, intent(in) :: nr,nri
36 real(8), intent(in) :: ri(nr),wcr(12,nr),rfmt(*)
37 integer, intent(in) :: ld
38 real(8), intent(out) :: grfmt(ld,3)
39 ! local variables
40 integer i
41 ! allocatable arrays
42 complex(8), allocatable :: zfmt(:),gzfmt(:,:)
43 allocate(zfmt(ld),gzfmt(ld,3))
44 ! convert real to complex spherical harmonic expansion
45 call rtozfmt(nr,nri,rfmt,zfmt)
46 ! compute the gradient
47 call gradzfmt(nr,nri,ri,wcr,zfmt,ld,gzfmt)
48 ! convert complex to real spherical harmonic expansion
49 do i=1,3
50  call ztorfmt(nr,nri,gzfmt(:,i),grfmt(:,i))
51 ! improve stability by smoothing the gradient
52  call rfmtsm(msmgmt,nr,nri,grfmt(:,i))
53 end do
54 deallocate(zfmt,gzfmt)
55 end subroutine
56 !EOC
57 
integer msmgmt
Definition: modmain.f90:222
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition: gradzfmt.f90:10
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition: rtozfmt.f90:7
subroutine rfmtsm(m, nr, nri, rfmt)
Definition: rfmtsm.f90:7
subroutine gradrfmt(nr, nri, ri, wcr, rfmt, ld, grfmt)
Definition: gradrfmt.f90:10
pure subroutine ztorfmt(nr, nri, zfmt, rfmt)
Definition: ztorfmt.f90:7