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