The Elk Code
genshtmat.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 General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: genshtmat
8 ! !INTERFACE:
9 subroutine genshtmat
10 ! !USES:
11 use modmain
12 ! !DESCRIPTION:
13 ! Generates the forward and backward spherical harmonic transformation (SHT)
14 ! matrices using the spherical covering set produced by the routine
15 ! {\tt sphcover}. These matrices are used to transform a function between its
16 ! $(l,m)$-expansion coefficients and its values at the $(\theta,\phi)$ points
17 ! on the sphere. Both real and complex SHT matrices are calculated and stored
18 ! in global arrays.
19 !
20 ! !REVISION HISTORY:
21 ! Created April 2003 (JKD)
22 !EOP
23 !BOC
24 implicit none
25 ! local variables
26 integer itp
27 real(8) v(3)
28 ! automatic arrays
29 real(8) tp(2,lmmaxo),vtp(3,lmmaxo),rlm(lmmaxo)
30 complex(8) ylm(lmmaxo)
31 !--------------------------------!
32 ! SHT matrices for lmaxi !
33 !--------------------------------!
34 ! allocate real SHT matrices
35 if (allocated(rbshti)) deallocate(rbshti)
36 allocate(rbshti(lmmaxi,lmmaxi))
37 if (allocated(rfshti)) deallocate(rfshti)
38 allocate(rfshti(lmmaxi,lmmaxi))
39 ! allocate complex SHT matrices
40 if (allocated(zbshti)) deallocate(zbshti)
41 allocate(zbshti(lmmaxi,lmmaxi))
42 if (allocated(zfshti)) deallocate(zfshti)
43 allocate(zfshti(lmmaxi,lmmaxi))
44 ! allocate single-precision complex copies
45 if (allocated(cbshti)) deallocate(cbshti)
46 allocate(cbshti(lmmaxi,lmmaxi))
47 if (allocated(cfshti)) deallocate(cfshti)
48 allocate(cfshti(lmmaxi,lmmaxi))
49 ! generate spherical covering set for lmaxi
50 call sphcover(lmmaxi,tp)
51 ! convert (theta, phi) angles to vectors
52 do itp=1,lmmaxi
53  call sctovec(tp(:,itp),vtp(:,itp))
54 end do
55 ! rotate the spherical covering set if required
56 if (trotsht) then
57  do itp=1,lmmaxi
58  v(1:3)=vtp(1:3,itp)
59  call r3mv(rotsht,v,vtp(:,itp))
60  end do
61 end if
62 ! generate real and complex spherical harmonics and set the backward SHT arrays
63 do itp=1,lmmaxi
64  call genrlmv(lmaxi,vtp(:,itp),rlm)
65  rbshti(itp,1:lmmaxi)=rlm(1:lmmaxi)
66  call genylmv(.false.,lmaxi,vtp(:,itp),ylm)
67  zbshti(itp,1:lmmaxi)=ylm(1:lmmaxi)
68 end do
69 ! find the forward SHT arrays
70 ! real
72 call rminv(lmmaxi,rfshti)
73 ! complex
75 call zminv(lmmaxi,zfshti)
76 ! make single-precision complex copies
79 !--------------------------------!
80 ! SHT matrices for lmaxo !
81 !--------------------------------!
82 ! allocate real SHT matrices
83 if (allocated(rbshto)) deallocate(rbshto)
84 allocate(rbshto(lmmaxo,lmmaxo))
85 if (allocated(rfshto)) deallocate(rfshto)
86 allocate(rfshto(lmmaxo,lmmaxo))
87 ! allocate complex SHT matrices
88 if (allocated(zbshto)) deallocate(zbshto)
89 allocate(zbshto(lmmaxo,lmmaxo))
90 if (allocated(zfshto)) deallocate(zfshto)
91 allocate(zfshto(lmmaxo,lmmaxo))
92 ! allocate single-precision complex copies
93 if (allocated(cbshto)) deallocate(cbshto)
94 allocate(cbshto(lmmaxo,lmmaxo))
95 if (allocated(cfshto)) deallocate(cfshto)
96 allocate(cfshto(lmmaxo,lmmaxo))
97 ! generate spherical covering set
98 call sphcover(lmmaxo,tp)
99 ! convert (theta, phi) angles to vectors
100 do itp=1,lmmaxo
101  call sctovec(tp(:,itp),vtp(:,itp))
102 end do
103 ! rotate the spherical covering set if required
104 if (trotsht) then
105  do itp=1,lmmaxo
106  v(1:3)=vtp(1:3,itp)
107  call r3mv(rotsht,v,vtp(:,itp))
108  end do
109 end if
110 ! generate real and complex spherical harmonics and set the backward SHT arrays
111 do itp=1,lmmaxo
112  call genrlmv(lmaxo,vtp(:,itp),rlm)
113  rbshto(itp,1:lmmaxo)=rlm(1:lmmaxo)
114  call genylmv(.false.,lmaxo,vtp(:,itp),ylm)
115  zbshto(itp,1:lmmaxo)=ylm(1:lmmaxo)
116 end do
117 ! find the forward SHT arrays
118 ! real
119 rfshto(1:lmmaxo,1:lmmaxo)=rbshto(1:lmmaxo,1:lmmaxo)
120 call rminv(lmmaxo,rfshto)
121 ! complex
122 zfshto(1:lmmaxo,1:lmmaxo)=zbshto(1:lmmaxo,1:lmmaxo)
123 call zminv(lmmaxo,zfshto)
124 ! make single-precision complex copies
125 cbshto(1:lmmaxo,1:lmmaxo)=zbshto(1:lmmaxo,1:lmmaxo)
126 cfshto(1:lmmaxo,1:lmmaxo)=zfshto(1:lmmaxo,1:lmmaxo)
127 
128 contains
129 
130 pure subroutine sctovec(tp,v)
131 implicit none
132 ! arguments
133 real(8), intent(in) :: tp(2)
134 real(8), intent(out) :: v(3)
135 ! local variables
136 real(8) t1
137 t1=sin(tp(1))
138 v(1)=t1*cos(tp(2))
139 v(2)=t1*sin(tp(2))
140 v(3)=cos(tp(1))
141 end subroutine
142 
143 end subroutine
144 !EOC
145 
subroutine genrlmv(lmax, v, rlm)
Definition: genrlmv.f90:10
complex(4), dimension(:,:), allocatable cfshti
Definition: modmain.f90:581
real(8), dimension(:,:), allocatable rbshto
Definition: modmain.f90:569
real(8), dimension(:,:), allocatable rbshti
Definition: modmain.f90:565
integer lmaxo
Definition: modmain.f90:201
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition: genylmv.f90:10
real(8), dimension(3, 3) rotsht
Definition: modmain.f90:563
subroutine genshtmat
Definition: genshtmat.f90:10
complex(8), dimension(:,:), allocatable zfshti
Definition: modmain.f90:575
integer lmmaxi
Definition: modmain.f90:207
real(8), dimension(:,:), allocatable rfshti
Definition: modmain.f90:567
logical trotsht
Definition: modmain.f90:561
complex(8), dimension(:,:), allocatable zfshto
Definition: modmain.f90:579
real(8), dimension(:,:), allocatable rfshto
Definition: modmain.f90:571
complex(4), dimension(:,:), allocatable cfshto
Definition: modmain.f90:582
complex(4), dimension(:,:), allocatable cbshti
Definition: modmain.f90:581
complex(8), dimension(:,:), allocatable zbshto
Definition: modmain.f90:577
subroutine zminv(n, a)
Definition: zminv.f90:7
subroutine rminv(n, a)
Definition: rminv.f90:7
complex(4), dimension(:,:), allocatable cbshto
Definition: modmain.f90:582
subroutine sphcover(n, tp)
Definition: sphcover.f90:10
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
pure subroutine sctovec(tp, v)
Definition: genshtmat.f90:131
complex(8), dimension(:,:), allocatable zbshti
Definition: modmain.f90:573
integer lmaxi
Definition: modmain.f90:205