-
Notifications
You must be signed in to change notification settings - Fork 103
/
SHReturnTapersM.f95
323 lines (278 loc) · 11.1 KB
/
SHReturnTapersM.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
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
subroutine SHReturnTapersM(theta0, lmax, m, tapers, eigenvalues, shannon, &
degrees, ntapers, exitstatus)
!------------------------------------------------------------------------------
!
! This subroutine will return all the eigenvalues and eigenfunctions for the
! space-concentration problem of a spherical cap of angular radius theta0.
! The returned eigenfunctions correspond to "geodesy" normalized spherical
! harmonic coefficients, and the eigenfunctions are further normalized such
! that they have unit power (i.e., the integral of the function squared over
! the sphere divided by 4 pi is 1, and the sum of the squares of their
! coefficients is 1). If the optional vector DEGREES is specified, then the
! eigenfunctions will be computed using only those degrees where DEGREES(l+1)
! is not zero.
!
! When possible, the eigenfunctions are calculated using the kernel of
! Grunbaum et al. 1982 and the eigenvalues are then calculated by integration
! using the definition of the space-concentration problem. Use of the
! Grunbaum et al. kernel is prefered over the space-concentration kernel as
! the eigenfunctions of the later are unreliable when there are several
! eigenvalues identical (within machine precision) to either 1 or zero. If,
! the optional parameter DEGREES is specified, and at least one element is
! zero for degrees greater or equal to abs(m), then the eigenfunctions and
! eigenvalues will instead be computed directly using the space-concentration
! kernel.
!
! Calling Parameters
!
! IN
! theta0 Angular radius of spherical cap in RADIANS.
! lmax Maximum spherical harmonic degree
! for the concentration problem.
! m Angular order of the concentration
! problem (m=0 corresponds to isotropic case).
!
! OUT
! tapers An (lmax+1) by (lmax+1-abs(m)) array containing
! all the eigenfunctions of the space-
! concentration kernel. Eigenfunctions
! are listed by columns in decreasing order
! corresponding to value of their eigenvalue.
! Only the first ntapers columns are non-zero.
! eigenvalues A vector of length lmax+1-abs(m) containing the
! eigenvalued corresponding to the individual
! eigenfunctions.
!
! OPTIONAL (IN)
! degrees Specify those degrees of the coupling matrix to
! compute. If degrees(l+1) is zero, degree l will not
! be computed.
!
! OPTIONAL (OUT)
! shannon Shannon number as calculated from the trace of the
! kernel.
! ntapers The number of non-zero tapers.
! 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.
!
! Copyright (c) 2005-2019, SHTOOLS
! All rights reserved.
!
!------------------------------------------------------------------------------
use SHTOOLS, only: ComputeDG82, ComputeDm, EigValVecSymTri, EigValVecSym, &
PreGLQ, PlmBar, PlmIndex
use ftypes
implicit none
real(dp), intent(in) :: theta0
integer(int32), intent(in) :: lmax, m
real(dp), intent(out) :: tapers(:,:), eigenvalues(:)
real(dp), intent(out), optional :: shannon
integer(int32), intent(in), optional :: degrees(:)
integer(int32), intent(out), optional :: ntapers
integer(int32), intent(out), optional :: exitstatus
integer(int32) :: l, n, n_int, j, i, astat(3), ind(lmax+1), use_dg82
real(dp) :: eval(lmax+1), pi, upper, lower, zero(lmax+1), w(lmax+1), h
real(dp), allocatable :: evec(:,:), dllmtri(:,:), p(:), dllm(:,:)
use_dg82 = 1
if (present(exitstatus)) exitstatus = 0
if (size(tapers(:,1)) < (lmax+1) .or. &
size(tapers(1,:)) < (lmax+1-abs(m)) ) then
print*, "Error --- SHReturnTapersM"
print*, "TAPERS must be dimensioned as (LMAX+1, LMAX+1-ABS(M)) " // &
"where LMAX is ", lmax
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(eigenvalues) < (lmax+1-abs(m)) ) then
print*, "Error --- SHReturnTapersM"
print*, "EIGENVALUES must be dimensioned as (LMAX+1-ABS(m)) " // &
"where LMAX is ", lmax
print*, "Input array is dimensioned as ", size(eigenvalues)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (m > lmax) then
print*, "Error --- SHReturnTapersM"
print*, "M must be less than or equal to LMAX."
print*, "M = ", m
print*, "LMAX = ", lmax
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
end if
if (present(degrees)) then
if (size(degrees) < lmax+1) then
print*, "Error --- SHReturnTapersM"
print*, "DEGREES must have dimension LMAX+1, where LMAX is ", lmax
print*, "Input array is dimensioned as ", size(degrees)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
end if
if (present(degrees)) then
if (sum(degrees(1+abs(m):lmax+1)) /= lmax + 1 - abs(m)) use_dg82 = 0
end if
allocate (evec(lmax+1, lmax+1-abs(m)), stat = astat(1))
allocate (dllmtri(lmax+1, lmax+1), stat = astat(2))
allocate (p((lmax+1)*(lmax+2)/2), stat = astat(3))
if (astat(1) /= 0 .or. astat(2) /= 0 .or. astat(3) /= 0) then
print*, "Error --- SHReturnTapersM"
print*, "Problem allocating arrays EVEC, DLLMTRI, and P", &
astat(1), astat(2), astat(3)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
pi = acos(-1.0_dp)
tapers = 0.0_dp
eigenvalues = 0.0_dp
eval = 0.0_dp
evec = 0.0_dp
if (use_dg82 == 1) then
!----------------------------------------------------------------------
!
! Calculate space-concentration kernel and the corresponding
! eigenfunctions of the Grunbaum et al. kernel.
!
!----------------------------------------------------------------------
n = lmax + 1 - abs(m)
if (present(exitstatus)) then
call ComputeDG82(dllmtri(1:n,1:n), lmax, m, theta0, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
call EigValVecSymTri(dllmtri(1:n,1:n), n, eval(1:n), &
tapers(1+abs(m):lmax+1,1:n), &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call ComputeDG82(dllmtri(1:n,1:n), lmax, m, theta0)
call EigValVecSymTri(dllmtri(1:n,1:n), n, eval(1:n), &
tapers(1+abs(m):lmax+1,1:n))
end if
!----------------------------------------------------------------------
!
! Calculate true eigenvalues
!
!----------------------------------------------------------------------
upper = 1.0_dp
lower = cos(theta0)
n_int = lmax + 1
if (present(exitstatus)) then
call PreGLQ(lower, upper, n_int, zero, w, exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call PreGLQ(lower, upper, n_int, zero, w)
end if
do i=1, n_int
if (present(exitstatus)) then
call PlmBar(p, lmax, zero(i), exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call PlmBar(p, lmax, zero(i))
end if
do j = 1, n
h = 0.0_dp
do l = abs(m), lmax
h = h + p(PlmIndex(l, abs(m))) * tapers(l+1, j)
end do
eigenvalues(j) = eigenvalues(j) + w(i) * h**2
end do
end do
if (m == 0) then
eigenvalues(1:n) = eigenvalues(1:n) / 2.0_dp
else
eigenvalues(1:n) = eigenvalues(1:n) / 4.0_dp
end if
else
!----------------------------------------------------------------------
!
! Calculate space-concentration kernel and the corresponding
! eigenfunctions.
!
!----------------------------------------------------------------------
ind = 0
i = 0
do l=abs(m), lmax, 1
if (degrees(l+1) /= 0) then
i = i + 1
ind(i) = l + 1
end if
end do
n = i
if (n /= 0) then
allocate (dllm(n, n), stat = astat(1))
if (astat(1) /= 0) then
print*, "Error --- SHReturnTapersM"
print*, "Problem allocating array DLLM ", astat(1)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
dllm = 0.0_dp
if (present(exitstatus)) then
call ComputeDm(dllmtri, lmax, m, theta0, degrees=degrees, &
exitstatus=exitstatus)
if (exitstatus /= 0) return
else
call ComputeDm(dllmtri, lmax, m, theta0, degrees=degrees)
end if
do i=1, n
do j=1, n
dllm(i, j) = dllmtri(ind(i), ind(j))
end do
end do
if (present(exitstatus)) then
call EigValVecSym(dllm(1:n, 1:n), n, &
eval(1:n), evec(1:n,1:n), &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call EigValVecSym(dllm, n, eval(1:n),&
evec(1:n,1:n))
end if
do i=1, n
tapers(ind(i), 1:n) = evec(i, 1:n)
end do
eigenvalues(1:n) = eval(1:n)
end if
end if
if (present(shannon)) then
shannon = sum(eigenvalues(1:lmax+1-abs(m)))
end if
if (present(ntapers)) then
ntapers = n
end if
! deallocate memory
call PlmBar (p, -1, zero(1))
deallocate (evec)
deallocate (dllmtri)
deallocate (p)
if (use_dg82 == 0 .and. n /= 0) deallocate (dllm)
end subroutine SHReturnTapersM