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