-
Notifications
You must be signed in to change notification settings - Fork 36
/
Copy pathclasses.py
2104 lines (1779 loc) · 74.2 KB
/
classes.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
# Copyright 2019 Pascal Audet & Helen Janiszewski
#
# This file is part of OBStools.
#
# Permission is hereby granted, free of charge, to any person obtaining a copy
# of this software and associated documentation files (the "Software"), to deal
# in the Software without restriction, including without limitation the rights
# to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
# copies of the Software, and to permit persons to whom the Software is
# furnished to do so, subject to the following conditions:
#
# The above copyright notice and this permission notice shall be included in
# all copies or substantial portions of the Software.
#
# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
# IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
# AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
# LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
# OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
# SOFTWARE.
import sys
from scipy.signal import stft, detrend
from scipy.linalg import norm
import matplotlib.pyplot as plt
import numpy as np
import pickle
from obspy.core import Stream, Trace, read
from obstools.atacr import utils, plotting
from pkg_resources import resource_filename
from pathlib import Path
np.seterr(all='ignore')
# np.set_printoptions(threshold=sys.maxsize)
class Power(object):
"""
Container for power spectra for each component, with any shape
Attributes
----------
c11 : :class:`~numpy.ndarray`
Power spectral density for component 1 (any shape)
c22 : :class:`~numpy.ndarray`
Power spectral density for component 2 (any shape)
cZZ : :class:`~numpy.ndarray`
Power spectral density for component Z (any shape)
cPP : :class:`~numpy.ndarray`
Power spectral density for component P (any shape)
"""
def __init__(self, c11=None, c22=None, cZZ=None, cPP=None):
self.c11 = c11
self.c22 = c22
self.cZZ = cZZ
self.cPP = cPP
class Cross(object):
"""
Container for cross-power spectra for each component pairs, with any shape
Attributes
----------
c12 : :class:`~numpy.ndarray`
Cross-power spectral density for components 1 and 2 (any shape)
c1Z : :class:`~numpy.ndarray`
Cross-power spectral density for components 1 and Z (any shape)
c1P : :class:`~numpy.ndarray`
Cross-power spectral density for components 1 and P (any shape)
c2Z : :class:`~numpy.ndarray`
Cross-power spectral density for components 2 and Z (any shape)
c2P : :class:`~numpy.ndarray`
Cross-power spectral density for components 2 and P (any shape)
cZP : :class:`~numpy.ndarray`
Cross-power spectral density for components Z and P (any shape)
"""
def __init__(self, c12=None, c1Z=None, c1P=None, c2Z=None, c2P=None,
cZP=None):
self.c12 = c12
self.c1Z = c1Z
self.c1P = c1P
self.c2Z = c2Z
self.c2P = c2P
self.cZP = cZP
class Rotation(object):
"""
Container for rotated spectra, with any shape
Attributes
----------
cHH : :class:`~numpy.ndarray`
Power spectral density for rotated horizontal component H (any shape)
cHZ : :class:`~numpy.ndarray`
Cross-power spectral density for components H and Z (any shape)
cHP : :class:`~numpy.ndarray`
Cross-power spectral density for components H and P (any shape)
coh : :class:`~numpy.ndarray`
Coherence between horizontal components
ph : :class:`~numpy.ndarray`
Phase of cross-power spectrum between horizontal components
tilt : float
Angle (azimuth) of tilt axis
coh_value : float
Maximum coherence
phase_value : float
Phase at maximum coherence
direc : :class:`~numpy.ndarray`
Directions for which the coherence is calculated
"""
def __init__(self, cHH=None, cHZ=None, cHP=None, coh=None, ph=None,
tilt=None, coh_value=None, phase_value=None, direc=None):
self.cHH = cHH
self.cHZ = cHZ
self.cHP = cHP
self.coh = coh
self.ph = ph
self.tilt = tilt
self.coh_value = coh_value
self.phase_value = phase_value
self.direc = direc
class DayNoise(object):
r"""
A DayNoise object contains attributes that associate
three-component raw (or deconvolved) traces, metadata information
and window parameters. The available methods carry out the quality
control steps and the average daily spectra for windows flagged as
"good".
Note
----
The object is initialized with :class:`~obspy.core.Trace` objects for
H1, H2, HZ and P components. Traces can be empty if data are not
available. Upon saving, those traces are discarded to save disk space.
Attributes
----------
tr1, tr2, trZ, trP : :class:`~obspy.core.Trace` object
Corresponding trace objects for components H1, H2, HZ and HP.
Traces can be empty (i.e., ``Trace()``) for missing components.
window : float
Length of time window in seconds
overlap : float
Fraction of overlap between adjacent windows
key : str
Station key for current object
dt : float
Sampling distance in seconds. Obtained from ``trZ`` object
npts : int
Number of points in time series. Obtained from ``trZ`` object
fs : float
Sampling frequency (in Hz). Obtained from ``trZ`` object
year : str
Year for current object (obtained from UTCDateTime). Obtained from
``trZ`` object
julday : str
Julian day for current object (obtained from UTCDateTime). Obtained
from ``trZ`` object
ncomp : int
Number of available components (either 2, 3 or 4). Obtained from
non-empty ``Trace`` objects
tf_list : Dict
Dictionary of possible transfer functions given the available
components.
Examples
--------
Get demo noise data as DayNoise object
>>> from obstools.atacr import DayNoise
>>> daynoise = DayNoise('demo')
Uploading demo data - March 04, 2012, station 7D.M08A
Now check its main attributes
>>> print(*[daynoise.tr1, daynoise.tr2, daynoise.trZ, daynoise.trP], sep="\n")
7D.M08A..BH1 | 2012-03-04T00:00:00.000000Z - 2012-03-04T23:59:59.800000Z | 5.0 Hz, 432000 samples
7D.M08A..BH2 | 2012-03-04T00:00:00.000000Z - 2012-03-04T23:59:59.800000Z | 5.0 Hz, 432000 samples
7D.M08A..BHZ | 2012-03-04T00:00:00.000000Z - 2012-03-04T23:59:59.800000Z | 5.0 Hz, 432000 samples
7D.M08A..BDH | 2012-03-04T00:00:00.000000Z - 2012-03-04T23:59:59.800000Z | 5.0 Hz, 432000 samples
>>> daynoise.window
7200.0
>>> daynoise.overlap
0.3
>>> daynoise.key
'7D.M08A'
>>> daynoise.ncomp
4
>>> daynoise.tf_list
{'ZP': True, 'Z1': True, 'Z2-1': True, 'ZP-21': True, 'ZH': True, 'ZP-H': True}
"""
def __init__(self, tr1=None, tr2=None, trZ=None, trP=None, window=7200.,
overlap=0.3, key=''):
# Load example data if initializing empty object
if tr1 == 'demo' or tr1 == 'Demo':
print("Uploading demo data - March 04, 2012, station 7D.M08A")
exmpl_path = Path(resource_filename('obstools', 'examples'))
fn = exmpl_path / 'data' / '2012.064*.SAC'
st = read(str(fn))
tr1 = st.select(component='1')[0]
tr2 = st.select(component='2')[0]
trZ = st.select(component='Z')[0]
trP = st.select(component='H')[0]
window = 7200.
overlap = 0.3
key = '7D.M08A'
# Check that all traces are valid Trace objects
for tr in [tr1, tr2, trZ, trP]:
if not isinstance(tr, Trace):
raise(Exception("Error initializing DayNoise object - "
+ str(tr)+" is not a Trace object"))
# Unpack everything
self.tr1 = tr1
self.tr2 = tr2
self.trZ = trZ
self.trP = trP
self.window = window
self.overlap = overlap
self.key = key
# Get trace attributes
zstats = self.trZ.stats
self.dt = zstats.delta
self.npts = zstats.npts
self.fs = zstats.sampling_rate
self.year = zstats.starttime.year
self.julday = zstats.starttime.julday
self.tkey = str(self.year) + '.' + str(self.julday)
# Get number of components for the available, non-empty traces
ncomp = np.sum(
[1 for tr in
Stream(traces=[tr1, tr2, trZ, trP]) if np.any(tr.data)])
self.ncomp = ncomp
# Build list of available transfer functions based on the number of
# components
if self.ncomp == 2:
self.tf_list = {'ZP': True, 'Z1': False, 'Z2-1': False,
'ZP-21': False, 'ZH': False, 'ZP-H': False}
elif self.ncomp == 3:
self.tf_list = {'ZP': False, 'Z1': True, 'Z2-1': True,
'ZP-21': False, 'ZH': True, 'ZP-H': False}
else:
self.tf_list = {'ZP': True, 'Z1': True, 'Z2-1': True,
'ZP-21': True, 'ZH': True, 'ZP-H': True}
self.QC = False
self.av = False
def QC_daily_spectra(self, pd=[0.004, 0.2], tol=1.5, alpha=0.05,
smooth=True, fig_QC=False, debug=False, save=None,
form='png'):
"""
Method to determine daily time windows for which the spectra are
anomalous and should be discarded in the calculation of the
transfer functions.
Parameters
----------
pd : list
Frequency corners of passband for calculating the spectra
tol : float
Tolerance threshold. If spectrum > std*tol, window is flagged as
bad
alpha : float
Confidence interval for f-test
smooth : boolean
Determines if the smoothed (True) or raw (False) spectra are used
fig_QC : boolean
Whether or not to produce a figure showing the results of the
quality control
debug : boolean
Whether or not to plot intermediate steps in the QC procedure
for debugging
save : :class:`~pathlib.Path` object
Relative path to figures folder
form : str
File format (e.g., 'png', 'jpg', 'eps')
Attributes
----------
ftX : :class:`~numpy.ndarray`
Windowed Fourier transform for the `X` component (can be either
1, 2, Z or P)
f : :class:`~numpy.ndarray`
Full frequency axis (Hz)
goodwins : list
List of booleans representing whether a window is good (True)
or not (False)
Examples
--------
Perform QC on DayNoise object using default values and plot final
figure
>>> from obstools.atacr import DayNoise
>>> daynoise = DayNoise('demo')
Uploading demo data - March 04, 2012, station 7D.M08A
>>> daynoise.QC_daily_spectra(fig_QC=True)
.. figure:: ../obstools/examples/figures/Figure_3a.png
:align: center
Print out new attribute of DayNoise object
>>> daynoise.goodwins
array([False, True, True, True, True, True, True, True, False,
False, True, True, True, True, True, True], dtype=bool)
"""
# Points in window
ws = int(self.window/self.dt)
# Number of points to overlap
ss = int(self.window*self.overlap/self.dt)
# hanning window
hanning = np.hanning(2*ss)
wind = np.ones(ws)
wind[0:ss] = hanning[0:ss]
wind[-ss:ws] = hanning[ss:ws]
# Get windowed Fourier transforms
ft1 = self.ft1 = None
ft2 = self.ft2 = None
ftZ = self.ftZ = None
ftP = self.ftP = None
# Calculate windowed FFTs and store as transpose
f, t, ftZ = stft(
self.trZ.data, self.fs, return_onesided=False, boundary=None,
padded=False, window=wind, nperseg=ws, noverlap=ss,
detrend='constant')
ftZ = ftZ * ws
self.ftZ = ftZ.T
if self.ncomp == 2 or self.ncomp == 4:
_f, _t, ftP = stft(
self.trP.data, self.fs, return_onesided=False, boundary=None,
padded=False, window=wind, nperseg=ws, noverlap=ss,
detrend='constant')
ftP = ftP * ws
self.ftP = ftP.T
if self.ncomp == 3 or self.ncomp == 4:
_f, _t, ft1 = stft(
self.tr1.data, self.fs, return_onesided=False, boundary=None,
padded=False, window=wind, nperseg=ws, noverlap=ss,
detrend='constant')
_f, _t, ft2 = stft(
self.tr2.data, self.fs, return_onesided=False, boundary=None,
padded=False, window=wind, nperseg=ws, noverlap=ss,
detrend='constant')
ft1 = ft1 * ws
ft2 = ft2 * ws
self.ft1 = ft1.T
self.ft2 = ft2.T
# Store frequency axis
self.f = f
# Get spectrograms for single day-long keys
psd1 = None
psd2 = None
psdZ = None
psdP = None
# Positive frequencies for PSD plots
faxis = int(len(f)/2)
f = f[0:faxis]
psdZ = np.abs(ftZ)**2*2./self.dt
psdZ = psdZ[0:faxis, :]
if self.ncomp == 2 or self.ncomp == 4:
psdP = np.abs(ftP)**2*2/self.dt
psdP = psdP[0:faxis, :]
if self.ncomp == 3 or self.ncomp == 4:
psd1 = np.abs(ft1)**2*2/self.dt
psd1 = psd1[0:faxis, :]
psd2 = np.abs(ft2)**2*2/self.dt
psd2 = psd2[0:faxis, :]
if fig_QC:
if self.ncomp == 2:
plt.figure(1)
plt.subplot(2, 1, 1)
plt.pcolormesh(t, f, np.log(psdZ), shading='auto')
plt.title('Z', fontdict={'fontsize': 8})
plt.subplot(2, 1, 2)
plt.pcolormesh(t, f, np.log(psdP), shading='auto')
plt.title('P', fontdict={'fontsize': 8})
plt.xlabel('Seconds')
plt.tight_layout()
if save:
fname = self.key + '.' + self.tkey + \
'.specgram_Z.P.' + form
if isinstance(save, Path):
fname = save / fname
plt.savefig(
str(fname), dpi=300, bbox_inches='tight', format=form)
else:
plt.show()
elif self.ncomp == 3:
plt.figure(1)
plt.subplot(3, 1, 1)
plt.pcolormesh(t, f, np.log(psd1), shading='auto')
plt.title('H1', fontdict={'fontsize': 8})
plt.subplot(3, 1, 2)
plt.pcolormesh(t, f, np.log(psd2), shading='auto')
plt.title('H2', fontdict={'fontsize': 8})
plt.subplot(3, 1, 3)
plt.pcolormesh(t, f, np.log(psdZ), shading='auto')
plt.title('Z', fontdict={'fontsize': 8})
plt.xlabel('Seconds')
plt.tight_layout()
if save:
fname = self.key + '.' + self.tkey + \
'.specgram_H1.H2.Z.' + form
if isinstance(save, Path):
fname = save / fname
plt.savefig(
str(fname), dpi=300, bbox_inches='tight', format=form)
else:
plt.show()
else:
plt.figure(1)
plt.subplot(4, 1, 1)
plt.pcolormesh(t, f, np.log(psd1), shading='auto')
plt.title('H1', fontdict={'fontsize': 8})
plt.subplot(4, 1, 2)
plt.pcolormesh(t, f, np.log(psd2), shading='auto')
plt.title('H2', fontdict={'fontsize': 8})
plt.subplot(4, 1, 3)
plt.pcolormesh(t, f, np.log(psdZ), shading='auto')
plt.title('Z', fontdict={'fontsize': 8})
plt.subplot(4, 1, 4)
plt.pcolormesh(t, f, np.log(psdP), shading='auto')
plt.title('P', fontdict={'fontsize': 8})
plt.xlabel('Seconds')
plt.tight_layout()
if save:
fname = self.key + '.' + self.tkey + \
'.specgram_H1.H2.Z.P.' + form
if isinstance(save, Path):
fname = save / fname
plt.savefig(
str(fname), dpi=300, bbox_inches='tight', format=form)
else:
plt.show()
# Select bandpass frequencies
ff = (f > pd[0]) & (f < pd[1])
if np.sum([(psd == 0.).any() for psd in
[psd1, psd2, psdZ, psdP] if psd is not None]) > 0.:
smooth = True
if smooth:
# Smooth out the log of the PSDs
sl_psd1 = None
sl_psd2 = None
sl_psdZ = None
sl_psdP = None
sl_psdZ = utils.smooth(np.log(psdZ, where=(psdZ > 0.)), 50, axis=0)
if self.ncomp == 2 or self.ncomp == 4:
sl_psdP = utils.smooth(
np.log(psdP, where=(psdP > 0.)), 50, axis=0)
if self.ncomp == 3 or self.ncomp == 4:
sl_psd1 = utils.smooth(
np.log(psd1, where=(psd1 > 0.)), 50, axis=0)
sl_psd2 = utils.smooth(
np.log(psd2, where=(psd2 > 0.)), 50, axis=0)
else:
# Take the log of the PSDs
sl_psd1 = None
sl_psd2 = None
sl_psdZ = None
sl_psdP = None
sl_psdZ = np.log(psdZ)
if self.ncomp == 2 or self.ncomp == 4:
sl_psdP = np.log(psdP)
if self.ncomp == 3 or self.ncomp == 4:
sl_psd1 = np.log(psd1)
sl_psd2 = np.log(psd2)
# Remove mean of the log PSDs
dsl_psdZ = sl_psdZ[ff, :] - np.mean(sl_psdZ[ff, :], axis=0)
if self.ncomp == 2:
dsl_psdP = sl_psdP[ff, :] - np.mean(sl_psdP[ff, :], axis=0)
dsls = [dsl_psdZ, dsl_psdP]
elif self.ncomp == 3:
dsl_psd1 = sl_psd1[ff, :] - np.mean(sl_psd1[ff, :], axis=0)
dsl_psd2 = sl_psd2[ff, :] - np.mean(sl_psd2[ff, :], axis=0)
dsls = [dsl_psd1, dsl_psd2, dsl_psdZ]
else:
dsl_psd1 = sl_psd1[ff, :] - np.mean(sl_psd1[ff, :], axis=0)
dsl_psd2 = sl_psd2[ff, :] - np.mean(sl_psd2[ff, :], axis=0)
dsl_psdP = sl_psdP[ff, :] - np.mean(sl_psdP[ff, :], axis=0)
dsls = [dsl_psd1, dsl_psd2, dsl_psdZ, dsl_psdP]
if self.ncomp == 2:
plt.figure(2)
plt.subplot(2, 1, 1)
plt.semilogx(f, sl_psdZ, 'g', lw=0.5)
plt.subplot(2, 1, 2)
plt.semilogx(f, sl_psdP, 'k', lw=0.5)
plt.tight_layout()
elif self.ncomp == 3:
plt.figure(2)
plt.subplot(3, 1, 1)
plt.semilogx(f, sl_psd1, 'r', lw=0.5)
plt.subplot(3, 1, 2)
plt.semilogx(f, sl_psd2, 'b', lw=0.5)
plt.subplot(3, 1, 3)
plt.semilogx(f, sl_psdZ, 'g', lw=0.5)
plt.tight_layout()
else:
plt.figure(2)
plt.subplot(4, 1, 1)
plt.semilogx(f, sl_psd1, 'r', lw=0.5)
plt.subplot(4, 1, 2)
plt.semilogx(f, sl_psd2, 'b', lw=0.5)
plt.subplot(4, 1, 3)
plt.semilogx(f, sl_psdZ, 'g', lw=0.5)
plt.subplot(4, 1, 4)
plt.semilogx(f, sl_psdP, 'k', lw=0.5)
plt.tight_layout()
if debug:
plt.show()
# Cycle through to kill high-std-norm windows
moveon = False
goodwins = np.repeat([True], len(t))
indwin = np.argwhere(goodwins == True)
while moveon == False:
ubernorm = np.empty((self.ncomp, np.sum(goodwins)))
for ind_u, dsl in enumerate(dsls):
normvar = np.zeros(np.sum(goodwins))
for ii, tmp in enumerate(indwin):
ind = np.copy(indwin)
ind = np.delete(ind, ii)
normvar[ii] = norm(np.std(dsl[:, ind], axis=1), ord=2)
ubernorm[ind_u, :] = np.median(normvar) - normvar
penalty = np.sum(ubernorm, axis=0)
plt.figure(4)
for i in range(self.ncomp):
plt.plot(range(0, np.sum(goodwins)), detrend(
ubernorm, type='constant')[i], 'o-')
if debug:
plt.show()
else:
plt.close('all')
plt.figure(5)
plt.plot(range(0, np.sum(goodwins)),
np.sum(ubernorm, axis=0), 'o-')
if debug:
plt.show()
else:
plt.close('all')
kill = penalty > tol*np.std(penalty)
if np.sum(kill) == 0:
self.goodwins = goodwins
moveon = True
trypenalty = penalty[np.argwhere(kill == False)].T[0]
if utils.ftest(penalty, 1, trypenalty, 1) < alpha:
goodwins[indwin[kill == True]] = False
indwin = np.argwhere(goodwins == True)
moveon = False
else:
moveon = True
self.goodwins = goodwins
if fig_QC:
power = Power(sl_psd1, sl_psd2, sl_psdZ, sl_psdP)
plot = plotting.fig_QC(
f, power, goodwins, self.ncomp, key=self.key)
# Save or show figure
if save:
fname = self.key + '.' + self.tkey + '.' + 'QC.' + form
if isinstance(save, Path):
fname = save / fname
plot.savefig(
str(fname), dpi=300, bbox_inches='tight', format=form)
else:
plot.show()
self.QC = True
def average_daily_spectra(self, calc_rotation=True, fig_average=False,
fig_coh_ph=False, save=None, form='png'):
"""
Method to average the daily spectra for good windows. By default, the
method will attempt to calculate the azimuth of maximum coherence
between horizontal components and the vertical component (for maximum
tilt direction), and use the rotated horizontals in the transfer
function calculations.
Parameters
----------
calc_rotation : boolean
Whether or not to calculate the tilt direction
fig_average : boolean
Whether or not to produce a figure showing the average daily
spectra
fig_coh_ph : boolean
Whether or not to produce a figure showing the maximum coherence
between H and Z
save : :class:`~pathlib.Path` object
Relative path to figures folder
form : str
File format (e.g., 'png', 'jpg', 'eps')
Attributes
----------
f : :class:`~numpy.ndarray`
Positive frequency axis for corresponding window parameters
power : :class:`~obstools.atacr.classes.Power`
Container for the Power spectra
cross : :class:`~obstools.atacr.classes.Cross`
Container for the Cross power spectra
rotation : :class:`~obstools.atacr.classes.Cross`, optional
Container for the Rotated power and cross spectra
Examples
--------
Average spectra for good windows using default values and plot final
figure
>>> from obstools.atacr import DayNoise
>>> daynoise = DayNoise('demo')
Uploading demo data - March 04, 2012, station 7D.M08A
>>> daynoise.QC_daily_spectra()
>>> daynoise.average_daily_spectra(fig_average=True)
.. figure:: ../obstools/examples/figures/Figure_3b.png
:align: center
Print out available attributes of DayNoise object
>>> daynoise.__dict__.keys()
dict_keys(['tr1', 'tr2', 'trZ', 'trP', 'window', 'overlap', 'key',
'dt', 'npts', 'fs', 'year', 'julday', 'tkey', 'ncomp', 'tf_list',
'QC', 'av', 'f', 'goodwins', 'power', 'cross', 'rotation'])
"""
if not self.QC:
print("Warning: processing daynoise object for " +
"QC_daily_spectra using default values")
self.QC_daily_spectra()
# Extract good windows
c11 = None
c22 = None
cZZ = None
cPP = None
cZZ = np.abs(
np.mean(
self.ftZ[self.goodwins, :]*np.conj(self.ftZ[self.goodwins, :]),
axis=0))
if self.ncomp == 2 or self.ncomp == 4:
cPP = np.abs(
np.mean(
self.ftP[self.goodwins, :]*np.conj(
self.ftP[self.goodwins, :]), axis=0))
if self.ncomp == 3 or self.ncomp == 4:
c11 = np.abs(
np.mean(
self.ft1[self.goodwins, :]*np.conj(
self.ft1[self.goodwins, :]), axis=0))
c22 = np.abs(
np.mean(
self.ft2[self.goodwins, :]*np.conj(
self.ft2[self.goodwins, :]), axis=0))
# Extract bad windows
bc11 = None
bc22 = None
bcZZ = None
bcPP = None
if np.sum(~self.goodwins) > 0:
bcZZ = np.abs(
np.mean(
self.ftZ[~self.goodwins, :]*np.conj(
self.ftZ[~self.goodwins, :]), axis=0))
if self.ncomp == 2 or self.ncomp == 4:
bcPP = np.abs(
np.mean(
self.ftP[~self.goodwins, :]*np.conj(
self.ftP[~self.goodwins, :]), axis=0))
if self.ncomp == 3 or self.ncomp == 4:
bc11 = np.abs(
np.mean(
self.ft1[~self.goodwins, :]*np.conj(
self.ft1[~self.goodwins, :]), axis=0))
bc22 = np.abs(
np.mean(
self.ft2[~self.goodwins, :]*np.conj(
self.ft2[~self.goodwins, :]), axis=0))
# Calculate mean of all good windows if component combinations exist
c12 = None
c1Z = None
c2Z = None
c1P = None
c2P = None
cZP = None
if self.ncomp == 3 or self.ncomp == 4:
c12 = np.mean(
self.ft1[self.goodwins, :] *
np.conj(self.ft2[self.goodwins, :]), axis=0)
c1Z = np.mean(
self.ft1[self.goodwins, :] *
np.conj(self.ftZ[self.goodwins, :]), axis=0)
c2Z = np.mean(
self.ft2[self.goodwins, :] *
np.conj(self.ftZ[self.goodwins, :]), axis=0)
if self.ncomp == 4:
c1P = np.mean(
self.ft1[self.goodwins, :] *
np.conj(self.ftP[self.goodwins, :]), axis=0)
c2P = np.mean(
self.ft2[self.goodwins, :] *
np.conj(self.ftP[self.goodwins, :]), axis=0)
if self.ncomp == 2 or self.ncomp == 4:
cZP = np.mean(
self.ftZ[self.goodwins, :] *
np.conj(self.ftP[self.goodwins, :]), axis=0)
# Store as attributes
self.power = Power(c11, c22, cZZ, cPP)
self.cross = Cross(c12, c1Z, c1P, c2Z, c2P, cZP)
bad = Power(bc11, bc22, bcZZ, bcPP)
if fig_average:
plot = plotting.fig_average(self.f, self.power, bad, self.goodwins,
self.ncomp, key=self.key)
if save:
fname = self.key + '.' + self.tkey + '.' + 'average.' + form
if isinstance(save, Path):
fname = save / fname
plot.savefig(
str(fname), dpi=300, bbox_inches='tight', format=form)
else:
plot.show()
if calc_rotation and self.ncomp >= 3:
cHH, cHZ, cHP, coh, ph, direc, tilt, coh_value, phase_value = \
utils.calculate_tilt(
self.ft1, self.ft2, self.ftZ, self.ftP, self.f,
self.goodwins)
self.rotation = Rotation(
cHH, cHZ, cHP, coh, ph, tilt, coh_value, phase_value, direc)
if fig_coh_ph:
plot = plotting.fig_coh_ph(coh, ph, direc)
# Save or show figure
if save:
fname = self.key + '.' + self.tkey + '.' + 'coh_ph.' + form
if isinstance(save, Path):
fname = save / fname
plot.savefig(
str(fname), dpi=300, bbox_inches='tight', format=form)
else:
plot.show()
else:
self.rotation = Rotation()
self.av = True
def save(self, filename):
"""
Method to save the object to file using `~Pickle`.
Parameters
----------
filename : str
File name
Examples
--------
Run demo through all methods
>>> from obstools.atacr import DayNoise
>>> daynoise = DayNoise('demo')
Uploading demo data - March 04, 2012, station 7D.M08A
>>> daynoise.QC_daily_spectra()
>>> daynoise.average_daily_spectra()
Save object
>>> daynoise.save('daynoise_demo.pkl')
Check that it has been saved
>>> import glob
>>> glob.glob("./daynoise_demo.pkl")
['./daynoise_demo.pkl']
"""
if not self.av:
print("Warning: saving before having calculated the average " +
"spectra")
# Remove original traces to save disk space
del self.tr1
del self.tr2
del self.trZ
del self.trP
file = open(str(filename), 'wb')
pickle.dump(self, file)
file.close()
class StaNoise(object):
"""
A StaNoise object contains attributes that associate
three-component raw (or deconvolved) traces, metadata information
and window parameters.
Note
----
The object is initially a container for
:class:`~obstools.atacr.classes.DayNoise` objects. Once the StaNoise
object is initialized (using the method `init()` or by calling the
`QC_sta_spectra()` method), each individual spectral quantity is unpacked
as an object attribute and the original `DayNoise` objects are removed
from memory. **DayNoise objects cannot be added or appended to the
StaNoise object once this is done**.
In addition, all spectral quantities associated with the
original `DayNoise` objects (now stored as attributes) are discarded as
the object is saved to disk and new container objects are defined and
saved.
Attributes
----------
daylist : list
A list of :class:`~obstools.atacr.classes.DayNoise` objects to process
and produce a station average
initialized : bool
Whether or not the object has been initialized - `False` unless one
of the methods have been called. When `True`, the `daylist` attribute
is deleted from memory
Examples
--------
Initialize empty object
>>> from obstools.atacr import StaNoise
>>> stanoise = StaNoise()
Initialize with DayNoise object
>>> from obstools.atacr import DayNoise
>>> daynoise = DayNoise('demo')
Uploading demo data - March 04, 2012, station 7D.M08A
>>> stanoise = StaNoise(daylist=[daynoise])
Add or append DayNoise object to StaNoise
>>> stanoise = StaNoise()
>>> stanoise += daynoise
>>> stanoise = StaNoise()
>>> stanoise.append(daynoise)
Import demo noise data with 4 DayNoise objects
>>> from obstools.atacr import StaNoise
>>> stanoise = StaNoise('demo')
Uploading demo data - March 01 to 04, 2012, station 7D.M08A
>>> stanoise.daylist
[<obstools.atacr.classes.DayNoise at 0x11e3ce8d0>,
<obstools.atacr.classes.DayNoise at 0x121c7ae10>,
<obstools.atacr.classes.DayNoise at 0x121ca5940>,
<obstools.atacr.classes.DayNoise at 0x121e7dd30>]
>>> stanoise.initialized
False
"""
def __init__(self, daylist=None):
def _load_dn(day):
exmpl_path = Path(resource_filename('obstools', 'examples'))
fn = '2012.'+day+'*.SAC'
fn = exmpl_path / 'data' / fn
st = read(str(fn))
tr1 = st.select(component='1')[0]
tr2 = st.select(component='2')[0]
trZ = st.select(component='Z')[0]
trP = st.select(component='H')[0]
window = 7200.
overlap = 0.3
key = '7D.M08A'
return DayNoise(tr1, tr2, trZ, trP, window, overlap, key)
self.daylist = []
self.initialized = False
self.QC = False
self.av = False
if isinstance(daylist, DayNoise):
daylist = [daylist]
elif daylist == 'demo' or daylist == 'Demo':
print("Uploading demo data - March 01 to 04, 2012, station " +
"7D.M08A")
self.daylist = [_load_dn('061'), _load_dn(
'062'), _load_dn('063'), _load_dn('064')]
if not daylist == 'demo' and daylist:
self.daylist.extend(daylist)
def __add__(self, other):
if isinstance(other, DayNoise):
other = StaNoise([other])
if not isinstance(other, StaNoise):
raise TypeError
daylist = self.daylist + other.daylist
return self.__class__(daylist=daylist)
def __iter__(self):
return List(self.daylist).__iter__()
def append(self, daynoise):
if isinstance(daynoise, DayNoise):
self.daylist.append(daynoise)
else:
msg = 'Append only supports a single DayNoise object as argument'
raise TypeError(msg)
return self
def extend(self, daynoise_list):
if isinstance(daynoise_list, list):
for _i in daynoise_list:
# Make sure each item in the list is a Grid object.
if not isinstance(_i, DayNoise):
msg = 'Extend only accepts a list of Daynoise objects.'
raise TypeError(msg)
self.daylist.extend(daynoise_list)
elif isinstance(daynoise_list, StaNoise):
self.daylist.extend(daynoise_list.daylist)
else:
msg = 'Extend only supports a list of DayNoise objects as ' +\
'argument.'
raise TypeError(msg)
return self
def init(self):
"""
Method to initialize the `StaNoise` object. This method is used to
unpack the spectral quantities from the original
:class:`~obstools.atacr.classes.DayNoise` objects and allow the
methods to proceed. The original
:class:`~obstools.atacr.classes.DayNoise` objects are deleted from
memory during this process.
Note
----
If the original :class:`~obstools.atacr.classes.DayNoise` objects