-
Notifications
You must be signed in to change notification settings - Fork 106
/
SHMultiTaperSE.f95
491 lines (417 loc) · 16.2 KB
/
SHMultiTaperSE.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
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
subroutine SHMultiTaperSE(mtse, sd, sh, lmax, tapers, taper_order, lmaxt, K, &
alpha, lat, lon, taper_wt, norm, csphase, exitstatus)
!------------------------------------------------------------------------------
!
! This subroutine will calculate the multitaper spectrum estimate utilizing
! the first K tapers. The standard error is calculated using an unbiased
! estimate of the sample variance.
!
! Calling Parameters
!
! IN
! sh Input spherical harmonic file.
! lmax Maximum degree of sh.
! tapers The eigenvector matrix returned from a program
! such as SHReturnTapers, where each column corresponds
! to the coefficients for a window for a single non-zero
! angular order.
! taper_order Angular order of the tapers.
! lmaxt Maximum degree of the eigentapers.
! K Number of tapers to use in the multitaper spectral
! estimation.
!
! OUT
! mtse Multitaper spectrum estimate, valid up to and including
! a maximum degree lmax-lmaxt.
! sd Standard error of the multitaper spectrum estimate.
!
! OPTIONAL (IN)
! alpha Euler angles used to rotate the localizing windows.
! lat Latitude to perform localized analysis (degrees).
! lon Longitude to perform localized analysis (degrees).
! taper_wt Weight to be applied to each direct spectral estimate.
! This should sum to unity.
! csphase: 1: Do not include the phase factor of (-1)^m
! -1: Apply the phase factor of (-1)^m.
! norm: Normalization to be used when calculating Legendre
! functions
! (1) "geodesy" (default)
! (2) Schmidt
! (3) unnormalized
! (4) orthonormalized
!
! 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.
!
! If the optional parameter alpha (or lat and lon) is input, then the
! spherical harmonic coefficients of the localizing windows will be
! rotated accordingly. To rotate a window originally centered at the
! north pole to (lat, long) given in degrees, use
!
! alpha(1) = 0.0
! alpha(2) = -(90.0_dp - lat)*pi/180.0_dp
! alpha(3) = -lon*pi/180.0_dp
!
! See documentation in file ShRotateRealCoef for further information on
! spherical harmonic rotations.
!
! Copyright (c) 2005-2019, SHTOOLS
! All rights reserved.
!
!------------------------------------------------------------------------------
use SHTOOLS, only: SHPowerSpectrum, SHRotateRealCoef, djpi2, &
CSPHASE_DEFAULT, SHGLQ, SHExpandGLQ, MakeGridGLQ
use ftypes
implicit none
real(dp), intent(out) :: mtse(:), sd(:)
real(dp), intent(in) :: sh(:,:,:), tapers(:,:)
integer, intent(in) :: lmax, lmaxt, K, taper_order(:)
real(dp), intent(in), optional :: alpha(:), lat, lon, taper_wt(:)
integer, intent(in), optional :: csphase, norm
integer, intent(out), optional :: exitstatus
integer :: i, l, phase, mnorm, astat(7), lmaxmul, nlat, nlong
real(dp) :: se(lmax-lmaxt+1, K), x(3), pi, factor
real(dp), allocatable, save :: zero(:), w(:)
integer, save :: first = 1, lmaxmul_last = -1
real(dp), allocatable :: shwin(:,:,:), shloc(:,:,:), dj(:,:,:), &
shwinrot(:,:,:), grid1glq(:,:), gridwinglq(:,:), &
temp(:,:)
!$OMP threadprivate(zero, w, first, lmaxmul_last)
if (present(exitstatus)) exitstatus = 0
pi = acos(-1.0_dp)
if (size(sh(:,1,1)) < 2 .or. size(sh(1,:,1)) < lmax+1 .or. &
size(sh(1,1,:)) < lmax+1) then
print*, "Error --- SHMultiTaperSE"
print*, "SH must be dimensioned (2,LMAX+1, LMAX+1) where LMAX is ", lmax
print*, "Input array is dimensioned ", size(sh(:,1,1)), &
size(sh(1,:,1)), size(sh(1,1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(tapers(:,1)) < lmaxt+1 .or. size(tapers(1,:)) < K) then
print*, "Error --- SHMultiTaperSE"
print*, "TAPERS must be dimensioned (LMAXT+1, K) where " // &
"LMAXT and K are, ", lmaxt, K
print*, "Input array is dimensioned ", size(tapers(:,1)), &
size(tapers(1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(taper_order) < K) then
print*, "Error --- SHMutltiTaperSE"
print*, "TAPER_ORDER must be dimensioned as (K) where K is ", K
print*, "Input dimension of array is ", size(taper_order)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(mtse) < lmax-lmaxt+1) then
print*, "Error --- SHMultiTaperSE"
print*, "MTSE must be dimensioned as (LMAX-LMAXT+1) where " // &
"LMAX and LMAXT are ", lmax, lmaxt
print*, "Input dimension of array is ", size(mtse)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(sd) < lmax-lmaxt+1) then
print*, "Error --- SHMultiTaperSE"
print*, "SD must be dimensioned as (LMAX-LMAXT+1) " // &
"where LMAX and LMAXT are ", lmax, lmaxt
print*, "Input dimension of array is ", size(sd)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (lmax < lmaxt) then
print*, "Error --- SHMultiTaperSE"
print*, "LMAX must be larger than LMAXT."
print*, "Input valuse of LMAX and LMAXT are ", lmax, lmaxt
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
end if
if (present(taper_wt)) then
if (size(taper_wt) < K) then
print*, "Error --- SHMultiTaperSE"
print*, "TAPER_WT must be dimensioned as (K) where K is ", K
print*, "Input dimension of array is ", size(taper_wt)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
end if
if (present(norm)) then
if (norm > 4 .or. norm < 1) then
print*, "Error --- SHMultiTaperSE"
print*, "Parameter NORM must be 1 (geodesy), 2 (Schmidt), " // &
"3 (unnormalized), or 4 (orthonormalized)."
print*, "Input value is ", norm
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
end if
mnorm = norm
else
mnorm = 1
end if
if (present(csphase)) then
if (csphase /= -1 .and. csphase /= 1) then
print*, "Error --- SHMultiTaperSE"
print*, "CSPHASE must be 1 (exclude) or -1 (include)."
print*, "Input value if ", csphase
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
else
phase = csphase
end if
else
phase = CSPHASE_DEFAULT
end if
if (present(lat) .and. .not. present(lon)) then
print*, "Error --- SHMultiTaperSE"
print*, "Both the optional parameters LAT and LON must be present."
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
end if
if (present(lon) .and. .not. present(lat)) then
print*, "Error --- SHMultiTaperSE"
print*, "Both the optional parameters LAT and LON must be present."
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
end if
if (present(alpha) .and. present(lat) .and. present(lon)) then
print*, "Error --- SHMultiTaperSE"
print*, "Only ALPHA or LAT and LON can be present, but not both"
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
end if
if (present(lat) .and. present(lon)) then
x(1) = 0.0_dp
x(2) = -(90.0_dp - lat) * pi / 180.0_dp
x(3) = -lon * pi / 180.0_dp
else if (present(alpha)) then
if (size(alpha) < 3) then
print*, "Error --- SHMultiTaperSE"
print*, "ALPHA must be dimensioned as (3)."
print*, "Input array is dimensioned as ", size(alpha)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
x = alpha
end if
lmaxmul = lmax + lmaxt
nlat = lmax + lmaxt + 1
nlong = 2 * (lmax + lmaxt) + 1
allocate (shwin(2,lmaxt+1,lmaxt+1), stat = astat(1))
allocate (shloc(2, lmax+lmaxt+1, lmax+lmaxt+1), stat = astat(2))
allocate (dj(lmaxt+1,lmaxt+1,lmaxt+1), stat = astat(3))
allocate (shwinrot(2,lmaxt+1,lmaxt+1), stat = astat(4))
allocate (grid1glq(nlat,nlong), stat = astat(5))
allocate (gridwinglq(nlat,nlong), stat = astat(6))
allocate (temp(nlat,nlong), stat = astat(7))
if (sum(astat(1:7)) /= 0) then
print*, "Error --- SHMultiTaperSE"
print*, "Problem allocating arrays SHWIN, SHLOC, DJ, " // &
"SHWINROT, GRID1GLQ, GRIDWINGLQ and TEMP", &
astat(1), astat(2), astat(3), astat(4), astat(5), &
astat(6), astat(7)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
if (first == 1) then
first = 0
lmaxmul_last = lmaxmul
allocate (zero(lmaxmul+1), stat = astat(1))
allocate (w(lmaxmul+1), stat = astat(2))
if (sum(astat(1:2)) /= 0) then
print*, "Error --- SHMultiTaperSE"
print*, "Problem allocating arrays ZERO and W", astat(1), astat(2)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
if (present(exitstatus)) then
call SHGLQ(lmaxmul, zero, w, csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHGLQ(lmaxmul, zero, w, csphase = phase, norm = mnorm)
end if
end if
if (lmaxmul /= lmaxmul_last) then
lmaxmul_last = lmaxmul
deallocate (zero)
deallocate (w)
allocate (zero(lmaxmul+1), stat = astat(1))
allocate (w(lmaxmul+1), stat = astat(2))
if (sum(astat(1:2)) /= 0) then
print*, "Error --- SHMultiTaperSE"
print*, "Problem allocating arrays ZERO and W", astat(1), astat(2)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
if (present(exitstatus)) then
call SHGLQ(lmaxmul, zero, w, csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHGLQ(lmaxmul, zero, w, csphase = phase, norm = mnorm)
end if
end if
mtse = 0.0_dp
sd = 0.0_dp
!--------------------------------------------------------------------------
!
! Calculate localized power spectra
!
!--------------------------------------------------------------------------
if (present(alpha) .or. (present(lat) .and. present(lon)) ) then
! Create rotation matrix used in the rotation routine.
if (present(exitstatus)) then
call djpi2(dj, lmaxt, exitstatus=exitstatus)
if (exitstatus /= 0) return
else
call djpi2(dj, lmaxt)
end if
end if
if (present(exitstatus)) then
call MakeGridGLQ(grid1glq, sh(1:2,1:lmax+1, 1:lmax+1), &
lmaxmul, zero = zero, csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call MakeGridGLQ(grid1glq, sh(1:2,1:lmax+1, 1:lmax+1), &
lmaxmul, zero = zero, csphase = phase, norm = mnorm)
end if
do i = 1, K
shwin = 0.0_dp
if (taper_order(i) < 0) then
shwin(2,1:lmaxt+1,abs(taper_order(i))+1) = tapers(1:lmaxt+1,i)
else
shwin(1,1:lmaxt+1,taper_order(i)+1) = tapers(1:lmaxt+1,i)
end if
if (present(alpha) .or. (present(lat) .and. present(lon)) ) then
if (present(exitstatus)) then
call SHRotateRealCoef(shwinrot, shwin, lmaxt, x, dj, &
exitstatus=exitstatus)
if (exitstatus /= 0) return
else
call SHRotateRealCoef(shwinrot, shwin, lmaxt, x, dj)
end if
shwin(1:2,1:lmaxt+1,1:lmaxt+1) = shwinrot(1:2,1:lmaxt+1,1:lmaxt+1)
end if
if (present(exitstatus)) then
call MakeGridGLQ(gridwinglq, shwin(1:2,1:lmaxt+1, 1:lmaxt+1), &
lmaxmul, zero = zero, csphase = phase, &
norm = mnorm, exitstatus = exitstatus)
if (exitstatus /= 0) return
temp(1:nlat,1:nlong) = grid1glq(1:nlat,1:nlong) &
* gridwinglq(1:nlat,1:nlong)
call SHExpandGLQ(shloc, lmaxmul, temp, w, zero = zero, &
csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
call SHPowerSpectrum(shloc, lmax-lmaxt, se(:,i), &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call MakeGridGLQ(gridwinglq, shwin(1:2,1:lmaxt+1, 1:lmaxt+1), &
lmaxmul, zero = zero, csphase = phase, &
norm = mnorm)
temp(1:nlat,1:nlong) = grid1glq(1:nlat,1:nlong) &
* gridwinglq(1:nlat,1:nlong)
call SHExpandGLQ(shloc, lmaxmul, temp, w, zero = zero, &
csphase = phase, norm = mnorm)
call SHPowerSpectrum(shloc, lmax-lmaxt, se(:,i))
end if
end do
if (present(taper_wt)) then
factor = sum(taper_wt(1:K))**2 - sum(taper_wt(1:K)**2)
factor = factor * sum(taper_wt(1:K))
factor = sum(taper_wt(1:K)**2) / factor
do l = 0, lmax-lmaxt, 1
mtse(l+1) = dot_product(se(l+1,1:K), taper_wt(1:K)) / &
sum(taper_wt(1:K))
if (K > 1) then
sd(l+1) = dot_product( (se(l+1,1:K) - mtse(l+1) )**2, &
taper_wt(1:K) ) * factor
end if
end do
else
do l = 0, lmax-lmaxt, 1
mtse(l+1) = sum(se(l+1,1:K)) / dble(K)
if (K > 1) then
sd(l+1) = sum( ( se(l+1,1:K) - mtse(l+1) )**2 ) / dble(K-1) &
/ dble(K) ! standard error!
end if
end do
end if
if (K > 1) sd(1:lmax-lmaxt+1) = sqrt(sd(1:lmax-lmaxt+1))
deallocate (shwin)
deallocate (shloc)
deallocate (dj)
deallocate (shwinrot)
deallocate (grid1glq)
deallocate (gridwinglq)
deallocate (temp)
end subroutine SHMultiTaperSE