The Elk Code
 
Loading...
Searching...
No Matches
rtozfmt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2013 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
6pure subroutine rtozfmt(nr,nri,rfmt,zfmt)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: nr,nri
11real(8), intent(in) :: rfmt(*)
12complex(8), intent(out) :: zfmt(*)
13! local variables
14integer i
15call rtozflmn(lmaxi,nri,lmmaxi,rfmt,zfmt)
16i=lmmaxi*nri+1
17call rtozflmn(lmaxo,nr-nri,lmmaxo,rfmt(i),zfmt(i))
18return
19
20contains
21
22!BOP
23! !ROUTINE: rtozflmn
24! !INTERFACE:
25pure subroutine rtozflmn(lmax,n,ld,rflm,zflm)
26! !INPUT/OUTPUT PARAMETERS:
27! lmax : maximum angular momentum (in,integer)
28! n : number of functions to convert (in,integer)
29! ld : leading dimension (in,integer)
30! rflm : coefficients of real spherical harmonic expansion (in,real(ld,n))
31! zflm : coefficients of complex spherical harmonic expansion
32! (out,complex(ld,n))
33! !DESCRIPTION:
34! Converts a real function, $r_{lm}$, expanded in terms of real spherical
35! harmonics into a complex spherical harmonic expansion, $z_{lm}$:
36! $$ z_{lm}=\begin{cases} \frac{1}{\sqrt{2}}(r_{lm}+i(-1)^mr_{l-m}) & m>0 \\
37! \frac{1}{\sqrt{2}}((-1)^mr_{l-m}-ir_{lm}) & m<0 \\
38! r_{lm} & m=0 \end{cases} $$
39! See routine {\tt genrlm}.
40!
41! !REVISION HISTORY:
42! Created April 2003 (JKD)
43!EOP
44!BOC
45implicit none
46! arguments
47integer, intent(in) :: lmax,n,ld
48real(8), intent(in) :: rflm(ld,n)
49complex(8), intent(out) :: zflm(ld,n)
50! local variables
51integer l,m,lm1,lm2
52! real constant 1/sqrt(2)
53real(8), parameter :: c1=0.7071067811865475244d0
54lm1=0
55do l=0,lmax
56 lm2=lm1+2*(l+1)
57 do m=-l,-1
58 lm1=lm1+1
59 lm2=lm2-1
60 if (mod(m,2) /= 0) then
61 zflm(lm1,1:n)=c1*cmplx(-rflm(lm2,1:n),-rflm(lm1,1:n),8)
62 else
63 zflm(lm1,1:n)=c1*cmplx(rflm(lm2,1:n),-rflm(lm1,1:n),8)
64 end if
65 end do
66 lm1=lm1+1
67 lm2=lm2-1
68 zflm(lm1,1:n)=rflm(lm1,1:n)
69 do m=1,l
70 lm1=lm1+1
71 lm2=lm2-1
72 if (mod(m,2) /= 0) then
73 zflm(lm1,1:n)=c1*cmplx(rflm(lm1,1:n),-rflm(lm2,1:n),8)
74 else
75 zflm(lm1,1:n)=c1*cmplx(rflm(lm1,1:n),rflm(lm2,1:n),8)
76 end if
77 end do
78end do
79end subroutine
80!EOC
81
82end subroutine
83
integer lmmaxi
Definition modmain.f90:207
integer lmaxo
Definition modmain.f90:201
integer lmaxi
Definition modmain.f90:205
integer lmmaxo
Definition modmain.f90:203
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition rtozfmt.f90:7
pure subroutine rtozflmn(lmax, n, ld, rflm, zflm)
Definition rtozfmt.f90:26