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