-
Notifications
You must be signed in to change notification settings - Fork 106
/
SHRotateTapers.f95
281 lines (246 loc) · 10.4 KB
/
SHRotateTapers.f95
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
subroutine SHRotateTapers(tapersrot, tapers, taper_order, lmax, nrot, x, dj, &
exitstatus)
!------------------------------------------------------------------------------
!
! This subroutine will rotate a set of real spherical cap tapers (originally
! centered at the North pole) corresponding to the angles listed in the
! input array X. Only the first NROT tapers are rotated, each of which is
! returned in a column of the output matrix TAPERSROT with spherical harmonic
! coefficients ordered according to SHCilmToVector.
!
! The rotation of a coordinate system or body can be viewed in two
! complementary ways involving three successive rotations. Both methods have
! the same initial and final configurations, and the angles listed in both
! schemes are the same.
!
! Scheme A: (I) Rotation about the z axis by alpha.
! (II) Rotation about the new y axis by beta.
! (III) Rotation about the new z axis by gamma.
!
! Scheme B: (I) Rotation about the z axis by gamma.
! (II) Rotation about the initial y axis by beta.
! (III) Rotation about the initial z axis by alpha.
!
! The rotations can further be viewed either as either a rotation of the
! coordinate system or the physical body.
!
! 1. Rotation of the coordinate system without rotation of the physical body,
! use x(alpha, beta, gamma).
!
! 2. Rotation of the physical body without rotation of the coordinate system,
! use x(-gamma, -beta, -alpha).
!
! To perform the inverse trasform of x(alpha, beta, gamma), use
! x(-gamma, -beta, -alpha).
!
! This routine uses the "y-convention" were rotations are about the y-axis
! instead of the x-axis.
!
! Calling Parameters
!
! IN
! tapers An (lmax+1) by nrot array containing
! all the eigenfunctions of the spherical cap
! concentration problem obtained from SHReturnTapers.
! The functions are listed by columns, ordered from
! best to worst concentrated.
! taper_order A vector of dimension (lmax+1)**2 denoting which
! order m corresponds to the column of tapers and
! eigenvalues.
! lmax Maximum spherical harmonic degree of the tapers.
! nrot Number of tapers to rotate.
! x Array or rotation angles in radians.
! dj Rotation matrix with dimension (lmax+1, lmax+1,
! lmax+1).
!
! OUT
! tapersrot Rotated real "geodesy" normalized spherical
! harmonic coefficients with dimension (lmax+1)**2
! by nrot.
!
! OPTIONAL (OUT)
! exitstatus If present, instead of executing a STOP when an error
! is encountered, the variable exitstatus will be
! returned describing the error.
! 0 = No errors;
! 1 = Improper dimensions of input array;
! 2 = Improper bounds for input variable;
! 3 = Error allocating memory;
! 4 = File IO error.
!
! Note: Before using this routine, please verify that the input Euler
! angles and signs give the expected results. Some people define the angle
! beta as a rotation with respect to the x axis.
!
! Copyright (c) 2005-2019, SHTOOLS
! All rights reserved.
!
!-------------------------------------------------------------------------------
use SHTOOLS, only: SHrtoc, SHctor, SHcilmtocindex, SHcindextocilm, &
SHRotateCoef, SHCilmToVector, CSPHASE_DEFAULT
use ftypes
implicit none
real(dp), intent(in) :: tapers(:,:), x(:), dj(:,:,:)
real(dp), intent(out) :: tapersrot(:,:)
integer, intent(in) :: taper_order(:), lmax, nrot
integer, intent(out), optional :: exitstatus
integer :: astat(5), i
real(dp), allocatable :: ccilm(:,:,:), cilm(:,:,:), cof(:,:), rcof(:,:), &
vec(:)
if (present(exitstatus)) exitstatus = 0
if (size(tapers(:,1)) < (lmax+1) .or. size(tapers(1,:)) < nrot) then
print*, "Error --- SHRotateTapers"
print*, "TAPERS must be dimensioned as ( LMAX+1, NROT ) " // &
"where LMAX = ", lmax, "and NROT = ", nrot
print*, "Input array is dimensioned as ", size(tapers(:,1)), &
size(tapers(1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(taper_order) < nrot) then
print*, "Error --- SHRotateTapers"
print*, "TAPER_ORDER must be dimensioned as NROT " // &
"where NROT = ", nrot
print*, "Input array is dimensioned as ", size(taper_order)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(tapersrot(:,1)) < (lmax+1)**2 .or. &
size(tapersrot(1,:)) < nrot) then
print*, "Error --- SHRotateTapers"
print*, "TAPERSROT must be dimensioned as ( (LMAX+1)**2, " // &
"NROT ), where LMAX = ", lmax, "and NROT = ", nrot
print*, "Input array is dimensioned as ", size(tapersrot(:,1)), &
size(tapersrot(1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(dj(:,1,1)) < lmax+1 .or. size(dj(1,:,1)) < lmax+1 &
.or. size(dj(1,1,:)) < lmax+1) then
print*, "Error --- SHRotateTapers"
print*, "DJ must be dimensioned as (LMAX+1, LMAX+1, LMAX+1) " // &
"where LMAX = ", lmax
print*, "Input array is dimensioned ", size(dj(:,1,1)), &
size(dj(1,:,1)), size(dj(1,1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(x) < 3) then
print*, "Error --- SHRotateTapers"
print*, "X must be dimensioned as (3)"
print*, "Input array is dimensioned ", size(x)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
allocate (ccilm(2,lmax+1,lmax+1), stat = astat(1))
allocate (cilm(2,lmax+1,lmax+1), stat = astat(2))
allocate (cof(2,((lmax+1)*(lmax+2))/2), stat = astat(3))
allocate (rcof(2,((lmax+1)*(lmax+2))/2), stat = astat(4))
allocate (vec((lmax+1)**2), stat = astat(5))
if (sum(astat(1:5)) /= 0) then
print*, "Error --- SHRotateTapers"
print*, "Problem allocating arrays CCILM, CILM, COF, RCOF and VEC", &
astat(1), astat(2), astat(3), astat(4), astat(5)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
do i=1, nrot, 1
cilm = 0.0_dp
if (taper_order(i) >= 0) then
cilm(1, 1:lmax+1, taper_order(i)+1) = tapers(1:lmax+1, i)
else
cilm(2, 1:lmax+1, abs(taper_order(i))+1) = tapers(1:lmax+1, i)
end if
if (CSPHASE_DEFAULT == 1) then
! Convert geodesy coefficients to Varshalovich et al. complex form
if (present(exitstatus)) then
call SHrtoc(cilm, ccilm, degmax=lmax, convention=2, &
switchcs=1, exitstatus=exitstatus)
if (exitstatus /= 0) return
else
call SHrtoc(cilm, ccilm, degmax=lmax, convention=2, switchcs=1)
end if
else
if (present(exitstatus)) then
call SHrtoc(cilm, ccilm, degmax=lmax, convention=2, &
switchcs=0, exitstatus=exitstatus)
if (exitstatus /= 0) return
else
call SHrtoc(cilm, ccilm, degmax=lmax, convention=2, switchcs=0)
end if
end if
if (present(exitstatus)) then
! Re-order complex coefficients to form a 2D vector
call SHcilmtocindex(ccilm, cof, lmax, exitstatus=exitstatus)
if (exitstatus /= 0) return
! Rotate complex re-ordered coefficients
call SHRotateCoef(x, cof, rcof, dj, lmax, exitstatus=exitstatus)
if (exitstatus /= 0) return
! Convert ordered complex coefficients back to a 3D array
call SHcindextocilm(rcof, ccilm, lmax, exitstatus=exitstatus)
if (exitstatus /= 0) return
else
! Re-order complex coefficients to form a 2D vector
call SHcilmtocindex(ccilm, cof, lmax)
! Rotate complex re-ordered coefficients
call SHRotateCoef(x, cof, rcof, dj, lmax)
! Convert ordered complex coefficients back to a 3D array
call SHcindextocilm(rcof, ccilm, lmax)
end if
if (CSPHASE_DEFAULT == 1) then
! Convert Varshalovich et al complex coefficients back to 4pi
! geodesy form
if (present(exitstatus)) then
call SHctor(ccilm, cilm, degmax=lmax, convention=2, &
switchcs=1, exitstatus=exitstatus)
if (exitstatus /= 0) return
else
call SHctor(ccilm, cilm, degmax=lmax, convention=2, &
switchcs=1)
end if
else
if (present(exitstatus)) then
call SHctor(ccilm, cilm, degmax=lmax, convention=2, &
switchcs=0, exitstatus=exitstatus)
if (exitstatus /= 0) return
else
call SHctor(ccilm, cilm, degmax=lmax, convention=2, &
switchcs=0)
end if
end if
if (present(exitstatus)) then
! Re-order real coefficients to form a 1-D vector
call SHCilmToVector(cilm, vec, lmax, exitstatus=exitstatus)
if (exitstatus /= 0) return
else
! Re-order real coefficients to form a 1-D vector
call SHCilmToVector(cilm, vec, lmax)
end if
tapersrot(1:(lmax+1)**2, i) = vec(1:(lmax+1)**2)
end do
deallocate (ccilm)
deallocate (cof)
deallocate (rcof)
deallocate (cilm)
deallocate (vec)
end subroutine SHRotateTapers