/
shwindow.py
1447 lines (1265 loc) · 58.2 KB
/
shwindow.py
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
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""
Spherical Harmonic Window classes
SHWindow: SHWindowCap, SHWindowMask
"""
from __future__ import absolute_import as _absolute_import
from __future__ import division as _division
from __future__ import print_function as _print_function
import numpy as _np
import matplotlib as _mpl
import matplotlib.pyplot as _plt
import copy as _copy
from .. import shtools as _shtools
from ..spectralanalysis import spectrum as _spectrum
from .shcoeffsgrid import SHCoeffs
from .shcoeffsgrid import SHGrid
__all__ = ['SHWindow', 'SHWindowCap', 'SHWindowMask']
class SHWindow(object):
"""
Class for spatio-spectral localization windows on the sphere.
The windows can be initialized from:
>>> x = SHWindow.from_cap(theta, lwin, [clat, clon, nwin])
>>> x = SHWindow.from_mask(SHGrid)
Each class instance defines the following class attributes:
kind : Either 'cap' or 'mask'.
tapers : Matrix containing the spherical harmonic coefficients
(in packed form) of either the unrotated spherical cap
localization windows or the localization windows
corresponding to the input mask.
coeffs : Array of spherical harmonic coefficients of the
rotated spherical cap localization windows. These are
'4pi' normalized and do not use the Condon-Shortley phase
factor.
eigenvalues : Concentration factors of the localization windows.
orders : The angular orders for each of the spherical cap
localization windows.
weights : Taper weights used with the multitaper spectral analyses.
Defaut is None.
lwin : Spherical harmonic bandwidth of the localization windows.
theta : Angular radius of the spherical cap localization domain
(default in degrees).
theta_degrees : True (default) if theta is in degrees.
nwin : Number of localization windows. Default is (lwin+1)**2.
nwinrot : The number of best concentrated windows that were rotated
and whose coefficients are stored in coeffs.
clat, clon : Latitude and longitude of the center of the rotated
spherical cap localization windows (default in degrees).
coord_degrees : True (default) if clat and clon are in degrees.
Each class instance provides the following methods:
to_array() : Return an array of the spherical harmonic
coefficients for taper i, where i=0 is the best
concentrated, optionally using a different
normalization convention.
to_shcoeffs() : Return the spherical harmonic coefficients of taper
i, where i=0 is the best concentrated, as a new
SHCoeffs class instance, optionally using a
different normalization convention.
to_shgrid() : Return as a new SHGrid instance a grid of taper i,
where i=0 is the best concentrated window.
number_concentrated() : Return the number of windows that have
concentration factors greater or equal to a
specified value.
degrees() : Return an array containing the spherical harmonic
degrees of the localization windows, from 0 to
lwin.
spectra() : Return the spectra of one or more localization
windows.
rotate() : Rotate the spherical cap tapers, originally located
at the north pole, to clat and clon and save the
spherical harmonic coefficients in the attribute
coeffs.
coupling_matrix() : Return the coupling matrix of the first nwin
localization windows.
biased_spectrum() : Calculate the multitaper (cross-) spectrum
expectation of a localized function.
multitaper_spectrum() : Return the multitaper power spectrum
estimate and uncertainty for the input
SHCoeffs class instance.
multitaper_cross_spectrum() : Return the multitaper cross-power
spectrum estimate and uncertainty for
two input SHCoeffs class instances.
copy() : Return a copy of the class instance.
plot_windows() : Plot the best concentrated localization windows
using a simple cylindrical projection.
plot_spectra() : Plot the spectra of the best-concentrated
localization windows.
plot_coupling_matrix() : Plot the multitaper coupling matrix.
info() : Print a summary of the data stored in the SHWindow
instance.
"""
def __init__(self):
"""Initialize with a factory method."""
pass
# ---- factory methods:
@classmethod
def from_cap(self, theta, lwin, clat=None, clon=None, nwin=None,
theta_degrees=True, coord_degrees=True, dj_matrix=None,
weights=None):
"""
Construct spherical cap localization windows.
Usage
-----
x = SHWindow.from_cap(theta, lwin, [clat, clon, nwin, theta_degrees,
coord_degrees, dj_matrix, weights])
Returns
-------
x : SHWindow class instance
Parameters
----------
theta : float
Angular radius of the spherical cap localization domain (default
in degrees).
lwin : int
Spherical harmonic bandwidth of the localization windows.
clat, clon : float, optional, default = None
Latitude and longitude of the center of the rotated spherical cap
localization windows (default in degrees).
nwin : int, optional, default (lwin+1)**2
Number of localization windows.
theta_degrees : bool, optional, default = True
True if theta is in degrees.
coord_degrees : bool, optional, default = True
True if clat and clon are in degrees.
dj_matrix : ndarray, optional, default = None
The djpi2 rotation matrix computed by a call to djpi2.
weights : ndarray, optional, default = None
Taper weights used with the multitaper spectral analyses.
"""
if theta_degrees:
tapers, eigenvalues, taper_order = _shtools.SHReturnTapers(
_np.radians(theta), lwin)
else:
tapers, eigenvalues, taper_order = _shtools.SHReturnTapers(
theta, lwin)
return SHWindowCap(theta, tapers, eigenvalues, taper_order,
clat, clon, nwin, theta_degrees, coord_degrees,
dj_matrix, weights, copy=False)
@classmethod
def from_mask(self, dh_mask, lwin, nwin=None, weights=None):
"""
Construct localization windows that are optimally concentrated within
the region specified by a mask.
Usage
-----
x = SHWindow.from_mask(dh_mask, lwin, [nwin, weights])
Returns
-------
x : SHWindow class instance
Parameters
----------
dh_mask :ndarray, shape (nlat, nlon)
A Driscoll and Healy (1994) sampled grid describing the
concentration region R. All elements should either be 1 (for inside
the concentration region) or 0 (for outside the concentration
region). The grid must have dimensions nlon=nlat or nlon=2*nlat,
where nlat is even.
lwin : int
The spherical harmonic bandwidth of the localization windows.
nwin : int, optional, default = (lwin+1)**2
The number of best concentrated eigenvalues and eigenfunctions to
return.
weights ndarray, optional, default = None
Taper weights used with the multitaper spectral analyses.
"""
if nwin is None:
nwin = (lwin + 1)**2
else:
if nwin > (lwin + 1)**2:
raise ValueError('nwin must be less than or equal to ' +
'(lwin + 1)**2. lwin = {:d} and nwin = {:d}'
.format(lwin, nwin))
if dh_mask.shape[0] % 2 != 0:
raise ValueError('The number of latitude bands in dh_mask ' +
'must be even. nlat = {:d}'
.format(dh_mask.shape[0]))
if (dh_mask.shape[1] != dh_mask.shape[0] and
dh_mask.shape[1] != 2 * dh_mask.shape[0]):
raise ValueError('dh_mask must be dimensioned as (n, n) or ' +
'(n, 2 * n). Input shape is ({:d}, {:d})'
.format(dh_mask.shape[0], dh_mask.shape[1]))
tapers, eigenvalues = _shtools.SHReturnTapersMap(dh_mask, lwin,
ntapers=nwin)
return SHWindowMask(tapers, eigenvalues, weights, copy=False)
def copy(self):
"""Return a deep copy of the class instance."""
return _copy.deepcopy(self)
def degrees(self):
"""
Return a numpy array listing the spherical harmonic degrees of the
localization windows from 0 to lwin.
Usage
-----
degrees = x.degrees()
Returns
-------
degrees : ndarray, shape (lwin+1)
numpy ndarray containing a list of the spherical harmonic degrees.
"""
return _np.arange(self.lwin + 1)
def number_concentrated(self, alpha):
"""
Return the number of localization windows that have concentration
factors greater or equal to alpha.
Usage
-----
k = x.number_concentrated(alpha)
Returns
-------
k : int
The number of windows with concentration factors greater or equal
to alpha.
Parameters
----------
alpha : float
The concentration factor, which is the power of the window within
the concentration region divided by the total power.
"""
return len(self.eigenvalues[self.eigenvalues >= alpha])
def to_array(self, itaper, normalization='4pi', csphase=1):
"""
Return the spherical harmonic coefficients of taper i as a numpy
array.
Usage
-----
coeffs = x.to_array(itaper, [normalization, csphase])
Returns
-------
coeffs : ndarray, shape (2, lwin+1, lwin+11)
3-D numpy ndarray of size of the spherical harmonic coefficients
of the window.
Parameters
----------
itaper : int
Taper number, where itaper=0 is the best concentrated.
normalization : str, optional, default = '4pi'
Normalization of the output coefficients: '4pi', 'ortho' or
'schmidt' for geodesy 4pi normalized, orthonormalized, or Schmidt
semi-normalized coefficients, respectively.
csphase : int, optional, default = 1
Condon-Shortley phase convention: 1 to exclude the phase factor,
or -1 to include it.
"""
if type(normalization) != str:
raise ValueError('normalization must be a string. ' +
'Input type was {:s}'
.format(str(type(normalization))))
if normalization.lower() not in ('4pi', 'ortho', 'schmidt'):
raise ValueError(
"normalization must be '4pi', 'ortho' " +
"or 'schmidt'. Provided value was {:s}"
.format(repr(output_normalization))
)
if csphase != 1 and csphase != -1:
raise ValueError(
"csphase must be 1 or -1. Input value was {:s}"
.format(repr(csphase))
)
return self._to_array(
itaper, normalization=normalization.lower(), csphase=csphase)
def to_shcoeffs(self, itaper, normalization='4pi', csphase=1):
"""
Return the spherical harmonic coefficients of taper i as a SHCoeffs
class instance.
Usage
-----
clm = x.to_shcoeffs(itaper, [normalization, csphase])
Returns
-------
clm : SHCoeffs class instance
Parameters
----------
itaper : int
Taper number, where itaper=0 is the best concentrated.
normalization : str, optional, default = '4pi'
Normalization of the output class: '4pi', 'ortho' or 'schmidt' for
geodesy 4pi-normalized, orthonormalized, or Schmidt semi-normalized
coefficients, respectively.
csphase : int, optional, default = 1
Condon-Shortley phase convention: 1 to exclude the phase factor,
or -1 to include it.
"""
if type(normalization) != str:
raise ValueError('normalization must be a string. ' +
'Input type was {:s}'
.format(str(type(normalization))))
if normalization.lower() not in set(['4pi', 'ortho', 'schmidt']):
raise ValueError(
"normalization must be '4pi', 'ortho' " +
"or 'schmidt'. Provided value was {:s}"
.format(repr(normalization))
)
if csphase != 1 and csphase != -1:
raise ValueError(
"csphase must be 1 or -1. Input value was {:s}"
.format(repr(csphase))
)
coeffs = self.to_array(itaper, normalization=normalization.lower(),
csphase=csphase)
return SHCoeffs.from_array(coeffs, normalization=normalization.lower(),
csphase=csphase, copy=False)
def to_shgrid(self, itaper, grid='DH2', zeros=None):
"""
Evaluate the coefficients of taper i on a spherical grid and return
a SHGrid class instance.
Usage
-----
f = x.to_shgrid(itaper, [grid, zeros])
Returns
-------
f : SHGrid class instance
Parameters
----------
itaper : int
Taper number, where itaper=0 is the best concentrated.
grid : str, optional, default = 'DH2'
'DH' or 'DH1' for an equisampled lat/lon grid with nlat=nlon, 'DH2'
for an equidistant lat/lon grid with nlon=2*nlat, or 'GLQ' for a
Gauss-Legendre quadrature grid.
zeros : ndarray, optional, default = None
The cos(colatitude) nodes used in the Gauss-Legendre Quadrature
grids.
Description
-----------
For more information concerning the spherical harmonic expansions and
the properties of the output grids, see the documentation for
SHExpandDH and SHExpandGLQ.
"""
if type(grid) != str:
raise ValueError('grid must be a string. ' +
'Input type was {:s}'
.format(str(type(grid))))
if grid.upper() in ('DH', 'DH1'):
gridout = _shtools.MakeGridDH(self.to_array(itaper), sampling=1,
norm=1, csphase=1)
return SHGrid.from_array(gridout, grid='DH', copy=False)
elif grid.upper() == 'DH2':
gridout = _shtools.MakeGridDH(self.to_array(itaper), sampling=2,
norm=1, csphase=1)
return SHGrid.from_array(gridout, grid='DH', copy=False)
elif grid.upper() == 'GLQ':
if zeros is None:
zeros, weights = _shtools.SHGLQ(self.lwin)
gridout = _shtools.MakeGridGLQ(self.to_array(itaper), zeros,
norm=1, csphase=1)
return SHGrid.from_array(gridout, grid='GLQ', copy=False)
else:
raise ValueError(
"grid must be 'DH', 'DH1', 'DH2', or 'GLQ'. " +
"Input value was {:s}".format(repr(grid)))
def multitaper_spectrum(self, clm, k, convention='power', unit='per_l',
**kwargs):
"""
Return the multitaper spectrum estimate and standard error.
Usage
-----
mtse, sd = x.multitaper_spectrum(clm, k, [convention, unit, lmax,
taper_wt, clat, clon,
coord_degrees])
Returns
-------
mtse : ndarray, shape (lmax-lwin+1)
The localized multitaper spectrum estimate, where lmax is the
spherical-harmonic bandwidth of clm, and lwin is the
spherical-harmonic bandwidth of the localization windows.
sd : ndarray, shape (lmax-lwin+1)
The standard error of the localized multitaper spectrum
estimate.
Parameters
----------
clm : SHCoeffs class instance
SHCoeffs class instance containing the spherical harmonic
coefficients of the global field to analyze.
k : int
The number of tapers to be utilized in performing the multitaper
spectral analysis.
convention : str, optional, default = 'power'
The type of output spectra: 'power' for power spectra, and
'energy' for energy spectra.
unit : str, optional, default = 'per_l'
The units of the output spectra. If 'per_l', the spectra contain
the total contribution for each spherical harmonic degree l. If
'per_lm', the spectra contain the average contribution for each
coefficient at spherical harmonic degree l.
lmax : int, optional, default = clm.lmax
The maximum spherical-harmonic degree of clm to use.
taper_wt : ndarray, optional, default = None
1-D numpy array of the weights used in calculating the multitaper
spectral estimates and standard error.
clat, clon : float, optional, default = 90., 0.
Latitude and longitude of the center of the spherical-cap
localization windows.
coord_degrees : bool, optional, default = True
True if clat and clon are in degrees.
"""
return self._multitaper_spectrum(clm, k, convention=convention,
unit=unit, **kwargs)
def multitaper_cross_spectrum(self, clm, slm, k, convention='power',
unit='per_l', **kwargs):
"""
Return the multitaper cross-spectrum estimate and standard error.
Usage
-----
mtse, sd = x.multitaper_cross_spectrum(clm, slm, k, [convention, unit,
lmax, taper_wt,
clat, clon,
coord_degrees])
Returns
-------
mtse : ndarray, shape (lmax-lwin+1)
The localized multitaper cross-spectrum estimate, where lmax is the
smaller of the two spherical-harmonic bandwidths of clm and slm,
and lwin is the spherical-harmonic bandwidth of the localization
windows.
sd : ndarray, shape (lmax-lwin+1)
The standard error of the localized multitaper cross-spectrum
estimate.
Parameters
----------
clm : SHCoeffs class instance
SHCoeffs class instance containing the spherical harmonic
coefficients of the first global field to analyze.
slm : SHCoeffs class instance
SHCoeffs class instance containing the spherical harmonic
coefficients of the second global field to analyze.
k : int
The number of tapers to be utilized in performing the multitaper
spectral analysis.
convention : str, optional, default = 'power'
The type of output spectra: 'power' for power spectra, and
'energy' for energy spectra.
unit : str, optional, default = 'per_l'
The units of the output spectra. If 'per_l', the spectra contain
the total contribution for each spherical harmonic degree l. If
'per_lm', the spectra contain the average contribution for each
coefficient at spherical harmonic degree l.
lmax : int, optional, default = min(clm.lmax, slm.lmax)
The maximum spherical-harmonic degree of the input coefficients
to use.
taper_wt : ndarray, optional, default = None
The weights used in calculating the multitaper cross-spectral
estimates and standard error.
clat, clon : float, optional, default = 90., 0.
Latitude and longitude of the center of the spherical-cap
localization windows.
coord_degrees : bool, optional, default = True
True if clat and clon are in degrees.
"""
return self._multitaper_cross_spectrum(clm, slm, k,
convention=convention,
unit=unit, **kwargs)
def biased_spectrum(self, power, k, convention='power', unit='per_l',
**kwargs):
"""
Calculate the multitaper (cross-)spectrum expectation of a
localized function.
Usage
-----
outspectrum = x.biased_spectrum(spectrum, k, [unit, power, taper_wt,
save_cg, ldata])
Returns
-------
outspectrum : ndarray, shape (ldata+lwin+1)
The expectation of the windowed spectrum, where ldata is the
spherical-harmonic bandwidth of the input spectrum, and lwin is the
spherical-harmonic bandwidth of the localization windows.
Parameters
----------
spectrum : ndarray, shape (ldata+1)
The global spectrum.
k : int
The number of best-concentrated localization windows to use in
constructing the windowed spectrum.
convention : str, optional, default = 'power'
The type of input global and output biased spectra: 'power' for
power spectra, and 'energy' for energy spectra.
unit : str, optional, default = 'per_l'
The units of the input global and output biased spectra. If
'per_l', the spectra contain the total contribution for each
spherical harmonic degree l. If 'per_lm', the spectra contain the
average contribution for each coefficient at spherical harmonic
degree l.
taper_wt : ndarray, optional, default = None
The weights used in calculating the multitaper spectral estimates
and standard error.
save_cg : int, optional, default = 0
If 1, the Clebsch-Gordon coefficients will be precomputed and saved
for future use. If 0, the Clebsch-Gordon coefficients will be
recomputed for each call.
ldata : int, optional, default = len(power)-1
The maximum degree of the global unwindowed spectrum.
"""
return self._biased_spectrum(power, k, convention=convention,
unit=unit, **kwargs)
def spectra(self, itaper=None, nwin=None, convention='power', unit='per_l',
base=10.):
"""
Return the spectra of one or more localization windows.
Usage
-----
spectra = x.spectra([itaper, nwin, convention, unit, base])
Returns
-------
spectra : ndarray, shape (lwin+1, nwin)
A matrix with each column containing the spectrum of a
localization window, and where the windows are arranged with
increasing concentration factors. If itaper is set, only a single
vector is returned, whereas if nwin is set, the first nwin spectra
are returned.
Parameters
----------
itaper : int, optional, default = None
The taper number of the output spectrum, where itaper=0
corresponds to the best concentrated taper.
nwin : int, optional, default = 1
The number of best concentrated localization window power spectra
to return.
convention : str, optional, default = 'power'
The type of spectrum to return: 'power' for power spectrum,
'energy' for energy spectrum, and 'l2norm' for the l2 norm
spectrum.
unit : str, optional, default = 'per_l'
If 'per_l', return the total contribution to the spectrum for each
spherical harmonic degree l. If 'per_lm', return the average
contribution to the spectrum for each coefficient at spherical
harmonic degree l. If 'per_dlogl', return the spectrum per log
interval dlog_a(l).
base : float, optional, default = 10.
The logarithm base when calculating the 'per_dlogl' spectrum.
Description
-----------
This function returns either the power spectrum, energy spectrum, or
l2-norm spectrum of one or more of the localization windows.
Total power is defined as the integral of the function squared over all
space, divided by the area the function spans. If the mean of the
function is zero, this is equivalent to the variance of the function.
The total energy is the integral of the function squared over all space
and is 4pi times the total power. The l2-norm is the sum of the
magnitude of the coefficients squared.
The output spectrum can be expresed using one of three units. 'per_l'
returns the contribution to the total spectrum from all angular orders
at degree l. 'per_lm' returns the average contribution to the total
spectrum from a single coefficient at degree l. The 'per_lm' spectrum
is equal to the 'per_l' spectrum divided by (2l+1). 'per_dlogl' returns
the contribution to the total spectrum from all angular orders over an
infinitessimal logarithmic degree band. The contrubution in the band
dlog_a(l) is spectrum(l, 'per_dlogl')*dlog_a(l), where a is the base,
and where spectrum(l, 'per_dlogl) is equal to
spectrum(l, 'per_l')*l*log(a).
"""
if itaper is None:
if nwin is None:
nwin = self.nwin
spectra = _np.zeros((self.lwin+1, nwin))
for iwin in range(nwin):
coeffs = self.to_array(iwin)
spectra[:, iwin] = _spectrum(coeffs, normalization='4pi',
convention=convention, unit=unit,
base=base)
else:
coeffs = self.to_array(itaper)
spectra = _spectrum(coeffs, normalization='4pi',
convention=convention, unit=unit, base=base)
return spectra
def coupling_matrix(self, lmax, nwin=None, weights=None, mode='full'):
"""
Return the coupling matrix of the first nwin tapers. This matrix
relates the global power spectrum to the expectation of the localized
multitaper spectrum.
Usage
-----
Mmt = x.coupling_matrix(lmax, [nwin, weights, mode])
Returns
-------
Mmt : ndarray, shape (lmax+lwin+1, lmax+1) or (lmax+1, lmax+1) or
(lmax-lwin+1, lmax+1)
Parameters
----------
lmax : int
Spherical harmonic bandwidth of the global power spectrum.
nwin : int, optional, default = x.nwin
Number of tapers used in the mutlitaper spectral analysis.
weights : ndarray, optional, default = x.weights
Taper weights used with the multitaper spectral analyses.
mode : str, opitonal, default = 'full'
'full' returns a biased output spectrum of size lmax+lwin+1. The
input spectrum is assumed to be zero for degrees l>lmax.
'same' returns a biased output spectrum with the same size
(lmax+1) as the input spectrum. The input spectrum is assumed to be
zero for degrees l>lmax.
'valid' returns a biased spectrum with size lmax-lwin+1. This
returns only that part of the biased spectrum that is not
influenced by the input spectrum beyond degree lmax.
"""
if weights is not None:
if nwin is not None:
if len(weights) != nwin:
raise ValueError(
'Length of weights must be equal to nwin. ' +
'len(weights) = {:d}, nwin = {:d}'.format(len(weights),
nwin))
else:
if len(weights) != self.nwin:
raise ValueError(
'Length of weights must be equal to nwin. ' +
'len(weights) = {:d}, nwin = {:d}'.format(len(weights),
self.nwin))
if mode == 'full':
return self._coupling_matrix(lmax, nwin=nwin, weights=weights)
elif mode == 'same':
cmatrix = self._coupling_matrix(lmax, nwin=nwin,
weights=weights)
return cmatrix[:lmax+1, :]
elif mode == 'valid':
cmatrix = self._coupling_matrix(lmax, nwin=nwin,
weights=weights)
return cmatrix[:lmax - self.lwin+1, :]
else:
raise ValueError("mode has to be 'full', 'same' or 'valid', not "
"{}".format(mode))
def plot_windows(self, nwin, show=True, fname=None):
"""
Plot the best-concentrated localization windows.
Usage
-----
x.plot_windows(nwin, [show, fname])
Parameters
----------
nwin : int
The number of localization windows to plot.
show : bool, optional, default = True
If True, plot the image to the screen.
fname : str, optional, default = None
If present, save the image to the file.
"""
if self.kind == 'cap':
if self.nwinrot is not None and self.nwinrot <= nwin:
nwin = self.nwinrot
maxcolumns = 5
ncolumns = min(maxcolumns, nwin)
nrows = _np.ceil(nwin / ncolumns).astype(int)
figsize = ncolumns * 2.4, nrows * 1.2 + 0.5
fig, axes = _plt.subplots(nrows, ncolumns, figsize=figsize)
for ax in axes[:-1, :].flatten():
for xlabel_i in ax.get_xticklabels():
xlabel_i.set_visible(False)
for ax in axes[:, 1:].flatten():
for ylabel_i in ax.get_yticklabels():
ylabel_i.set_visible(False)
for itaper in range(min(self.nwin, nwin)):
evalue = self.eigenvalues[itaper]
ax = axes.flatten()[itaper]
gridout = _shtools.MakeGridDH(self.to_array(itaper), sampling=2,
norm=1, csphase=1)
ax.imshow(gridout, origin='upper',
extent=(0., 360., -90., 90.))
ax.set_title('concentration: {:2.2f}'.format(evalue))
fig.tight_layout(pad=0.5)
if show:
_plt.show()
if fname is not None:
fig.savefig(fname)
return fig, axes
def plot_spectra(self, nwin, convention='power', unit='per_l', base=10.,
xscale='lin', yscale='log', show=True, fname=None):
"""
Plot the spectra of the best-concentrated localization windows.
Usage
-----
x.plot_spectra(nwin, [convention, unit, base, xscale, yscale,
show, fname])
Parameters
----------
nwin : int
The number of localization windows to plot.
convention : str, optional, default = 'power'
The type of spectra to plot: 'power' for power spectrum, and
'energy' for energy spectrum.
unit : str, optional, default = 'per_l'
If 'per_l', return the total contribution to the spectrum for each
spherical harmonic degree l. If 'per_lm', return the average
contribution to the spectrum for each coefficient at spherical
harmonic degree l. If 'per_dlogl', return the spectrum per log
interval dlog_a(l).
base : float, optional, default = 10.
The logarithm base when calculating the 'per_dlogl' spectrum.
xscale : str, optional, default = 'lin'
Scale of the x axis: 'lin' for linear or 'log' for logarithmic.
yscale : str, optional, default = 'log'
Scale of the y axis: 'lin' for linear or 'log' for logarithmic.
show : bool, optional, default = True
If True, plot the image to the screen.
fname : str, optional, default = None
If present, save the image to the file.
"""
degrees = self.degrees()
spectrum = self.spectra(nwin=nwin, convention=convention, unit=unit,
base=base)
maxcolumns = 5
ncolumns = min(maxcolumns, nwin)
nrows = _np.ceil(nwin / ncolumns).astype(int)
figsize = ncolumns * 2.4, nrows * 1.2 + 0.5
fig, axes = _plt.subplots(nrows, ncolumns, figsize=figsize)
for ax in axes[:-1, :].flatten():
for xlabel_i in ax.get_xticklabels():
xlabel_i.set_visible(False)
for ax in axes[:, 1:].flatten():
for ylabel_i in ax.get_yticklabels():
ylabel_i.set_visible(False)
for itaper in range(min(self.nwin, nwin)):
evalue = self.eigenvalues[itaper]
ax = axes.flatten()[itaper]
ax.set_xlabel('degree l')
if (convention == 'power'):
ax.set_ylabel('power')
else:
ax.set_ylabel('energy')
if yscale == 'log':
ax.set_yscale('log', basey=base)
if xscale == 'log':
ax.set_xscale('log', basex=base)
ax.plot(degrees[1:], spectrum[1:, itaper])
else:
ax.plot(degrees[0:], spectrum[0:, itaper])
ax.grid(True, which='both')
ax.set_title('concentration: {:2.2f}'.format(evalue))
fig.tight_layout(pad=0.5)
if show:
_plt.show()
if fname is not None:
fig.savefig(fname)
return fig, axes
def plot_coupling_matrix(self, lmax, nwin=None, weights=None, mode='full',
show=True, fname=None):
"""
Plot the multitaper coupling matrix.
This matrix relates the global power spectrum to the expectation of
the localized multitaper spectrum.
Usage
-----
x.plot_coupling_matrix(lmax, [nwin, weights, mode, show, fname])
Parameters
----------
lmax : int
Spherical harmonic bandwidth of the global power spectrum.
nwin : int, optional, default = x.nwin
Number of tapers used in the mutlitaper spectral analysis.
weights : ndarray, optional, default = x.weights
Taper weights used with the multitaper spectral analyses.
mode : str, opitonal, default = 'full'
'full' returns a biased output spectrum of size lmax+lwin+1. The
input spectrum is assumed to be zero for degrees l>lmax.
'same' returns a biased output spectrum with the same size
(lmax+1) as the input spectrum. The input spectrum is assumed to be
zero for degrees l>lmax.
'valid' returns a biased spectrum with size lmax-lwin+1. This
returns only that part of the biased spectrum that is not
influenced by the input spectrum beyond degree lmax.
show : bool, optional, default = True
If True, plot the image to the screen.
fname : str, optional, default = None
If present, save the image to the file.
"""
figsize = _mpl.rcParams['figure.figsize']
figsize[0] = figsize[1]
fig = _plt.figure(figsize=figsize)
ax = fig.add_subplot(111)
ax.imshow(self.coupling_matrix(lmax, nwin=nwin, weights=weights,
mode=mode), aspect='auto')
ax.set_xlabel('input power') # matrix index 1 (columns)
ax.set_ylabel('output power') # matrix index 0 (rows)
fig.tight_layout(pad=0.1)
if show:
_plt.show()
if fname is not None:
fig.savefig(fname)
return fig, ax
def info(self):
"""
Print a summary of the data stored in the SHWindow instance.
Usage
-----
x.info()
"""
self._info()
class SHWindowCap(SHWindow):
"""Class for localization windows concentrated within a spherical cap."""
@staticmethod
def istype(kind):
return kind == 'cap'
def __init__(self, theta, tapers, eigenvalues, taper_order,
clat, clon, nwin, theta_degrees, coord_degrees, dj_matrix,
weights, copy=True):
self.kind = 'cap'
self.theta = theta
self.clat = clat
self.clon = clon
self.lwin = tapers.shape[0] - 1
self.theta_degrees = theta_degrees
self.coord_degrees = coord_degrees
self.dj_matrix = dj_matrix
self.weights = weights
self.nwinrot = None
if nwin is not None:
self.nwin = nwin
else:
self.nwin = tapers.shape[1]
if self.nwin > (self.lwin + 1)**2:
raise ValueError('nwin must be less than or equal to ' +
'(lwin+1)**2. nwin = {:s} and lwin = {:s}.'
.format(repr(self.nwin), repr(self.lwin)))
if copy:
self.tapers = _np.copy(tapers[:, :self.nwin])
self.eigenvalues = _np.copy(eigenvalues[:self.nwin])
self.orders = _np.copy(taper_order[:self.nwin])
else:
self.tapers = tapers[:, :self.nwin]
self.eigenvalues = eigenvalues[:self.nwin]
self.orders = taper_order[:self.nwin]
# If the windows aren't rotated, don't store them.
if self.clat is None and self.clon is None:
self.coeffs = None
else:
self.rotate(clat=self.clat, clon=self.clon,
coord_degrees=self.coord_degrees,
dj_matrix=self.dj_matrix)
def _taper2coeffs(self, itaper):
"""
Return the spherical harmonic coefficients of the unrotated taper i
as an array, where i = 0 is the best concentrated.
"""
taperm = self.orders[itaper]
coeffs = _np.zeros((2, self.lwin + 1, self.lwin + 1))
if taperm < 0:
coeffs[1, :, abs(taperm)] = self.tapers[:, itaper]
else:
coeffs[0, :, abs(taperm)] = self.tapers[:, itaper]
return coeffs
def _to_array(self, itaper, normalization='4pi', csphase=1):
"""
Return the spherical harmonic coefficients of taper i as an
array, where i = 0 is the best concentrated.
"""
if self.coeffs is None:
coeffs = _np.copy(self._taper2coeffs(itaper))
else:
if itaper > self.nwinrot - 1:
raise ValueError('itaper must be less than or equal to ' +
'nwinrot - 1. itaper = {:d}, nwinrot = {:d}'
.format(itaper, self.nwinrot))
coeffs = _shtools.SHVectorToCilm(self.coeffs[:, itaper])
if normalization == 'schmidt':
for l in range(self.lwin + 1):
coeffs[:, l, :l+1] *= _np.sqrt(2.0 * l + 1.0)
elif normalization == 'ortho':
coeffs *= _np.sqrt(4.0 * _np.pi)
if csphase == -1:
for m in range(self.lwin + 1):
if m % 2 == 1:
coeffs[:, :, m] = - coeffs[:, :, m]
return coeffs
def rotate(self, clat, clon, coord_degrees=True, dj_matrix=None,
nwinrot=None):
""""
Rotate the spherical-cap windows centered on the north pole to clat
and clon, and save the spherical harmonic coefficients in the
attribute coeffs.
Usage
-----
x.rotate(clat, clon [coord_degrees, dj_matrix, nwinrot])
Parameters
----------
clat, clon : float
Latitude and longitude of the center of the rotated spherical-cap
localization windows (default in degrees).
coord_degrees : bool, optional, default = True
True if clat and clon are in degrees.
dj_matrix : ndarray, optional, default = None
The djpi2 rotation matrix computed by a call to djpi2.
nwinrot : int, optional, default = (lwin+1)**2