The Elk Code
genrlmv.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: genrlmv
8 ! !INTERFACE:
9 subroutine genrlmv(lmax,v,rlm)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! lmax : maximum angular momentum (in,integer)
12 ! v : input vector (in,real(3))
13 ! rlm : array of real spherical harmonics (out,real((lmax+1)**2))
14 ! !DESCRIPTION:
15 ! Generates a sequence of real spherical harmonics evaluated at angles
16 ! $(\theta,\phi)$ for $0<l<l_{\rm max}$. The values are returned in a packed
17 ! array {\tt rlm} indexed with $j=l(l+1)+m+1$. Real spherical harmonics are
18 ! defined by
19 ! $$ R_{lm}(\theta,\phi)=\begin{cases}
20 ! \sqrt{2}\,{\rm Re}\{Y_{lm}(\theta,\phi)\} & m>0 \\
21 ! \sqrt{2}\,{\rm Im}\{Y_{lm}(\theta,\phi)\} & m<0 \\
22 ! {\rm Re}\{Y_{lm}(\theta,\phi)\} & m=0
23 ! \end{cases} $$
24 ! where $Y_{lm}$ are the complex spherical harmonics. These functions are
25 ! orthonormal and complete and may be used for expanding real-valued functions
26 ! on the sphere. This routine is numerically stable and accurate to near
27 ! machine precision for $l\le 50$. See routine {\tt genylmv}.
28 !
29 ! !REVISION HISTORY:
30 ! Created March 2004 (JKD)
31 !EOP
32 !BOC
33 implicit none
34 ! arguments
35 integer, intent(in) :: lmax
36 real(8), intent(in) :: v(3)
37 real(8), intent(out) :: rlm(*)
38 ! local variables
39 integer l,lm,n
40 real(8), parameter :: sqtwo=1.4142135623730950488d0
41 ! automatic arrays
42 complex(8) ylm((lmax+1)**2)
43 if ((lmax < 0).or.(lmax > 50)) then
44  write(*,*)
45  write(*,'("Error(genrlmv): lmax out of range : ",I8)') lmax
46  write(*,*)
47  stop
48 end if
49 ! generate complex spherical harmonics
50 call genylmv(.false.,lmax,v,ylm)
51 ! convert to real spherical harmonics
52 rlm(1)=dble(ylm(1))
53 do l=1,lmax
54  n=l-1
55  lm=l**2+1
56  rlm(lm:lm+n)=sqtwo*aimag(ylm(lm:lm+n))
57  lm=lm+l
58  rlm(lm)=dble(ylm(lm))
59  lm=lm+1
60  rlm(lm:lm+n)=sqtwo*dble(ylm(lm:lm+n))
61 end do
62 end subroutine
63 !EOC
64 
subroutine genrlmv(lmax, v, rlm)
Definition: genrlmv.f90:10
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition: genylmv.f90:10