-
Notifications
You must be signed in to change notification settings - Fork 24
/
Copy pathbenchmark_cap_vs_mtuq.py
235 lines (187 loc) · 8.33 KB
/
benchmark_cap_vs_mtuq.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
import os
import numpy as np
from mtuq import read, open_db, download_greens_tensors
from mtuq.event import MomentTensor
from mtuq.graphics import plot_waveforms2, plot_beachball, plot_misfit_dc
from mtuq.grid import DoubleCoupleGridRegular
from mtuq.grid_search import grid_search
from mtuq.misfit import Misfit
from mtuq.process_data import ProcessData
from mtuq.util import fullpath, merge_dicts, save_json
from mtuq.util.cap import parse_station_codes, Trapezoid
if __name__=='__main__':
#
# Given seven "fundamental" moment tensors, generates MTUQ synthetics and
# compares with corresponding CAP/FK synthetics
#
# Before running this script, it is necessary to unpack the CAP/FK
# synthetics using data/tests/unpack.bash
#
# This script is similar to examples/SerialGridSearch.DoubleCouple.py,
# except here we consider only seven grid points rather than an entire
# grid, and here the final plots are a comparison of MTUQ and CAP/FK
# synthetics rather than a comparison of data and synthetics
#
# Because of the idiosyncratic way CAP implements source-time function
# convolution, it's not expected that CAP and MTUQ synthetics will match
# exactly. CAP's "conv" function results in systematic magnitude-
# dependent shifts between origin times and arrival times. We deal with
# this by applying magnitude-dependent time-shifts to MTUQ synthetics
# (which normally lack such shifts) at the end of the benchmark. Even with
# this correction, the match will not be exact because CAP applies the
# shifts before tapering and MTUQ after tapering. The resulting mismatch
# will usually be apparent in body-wave windows, but not in surface-wave
# windows
#
# Note that CAP works with dyne,cm and MTUQ works with N,m, so to make
# comparisons we convert CAP output from the former to the latter
#
# The CAP/FK synthetics used in the comparison were generated by
# uafseismo/capuaf:46dd46bdc06e1336c3c4ccf4f99368fe99019c88
# using the following commands
#
# source #0 (explosion):
# cap.pl -H0.02 -P1/15/60 -p1 -S2/10/0 -T15/150 -D1/1/0.5 -C0.1/0.333/0.025/0.0625 -Y1 -Zweight_test.dat -Mscak_34 -m4.5 -I1 -R0/1.178/90/45/90 20090407201255351
#
# source #1 (on-diagonal)
# cap.pl -H0.02 -P1/15/60 -p1 -S2/10/0 -T15/150 -D1/1/0.5 -C0.1/0.333/0.025/0.0625 -Y1 -Zweight_test.dat -Mscak_34 -m4.5 -I1 -R-0.333/0.972/90/45/90 20090407201255351
#
# source #2 (on-diagonal)
# cap.pl -H0.02 -P1/15/60 -p1 -S2/10/0 -T15/150 -D1/1/0.5 -C0.1/0.333/0.025/0.0625 -Y1 -Zweight_test.dat -Mscak_34 -m4.5 -I1 -R-0.333/0.972/45/90/180 20090407201255351
#
# source #3 (on-diagonal):
# cap.pl -H0.02 -P1/15/60 -p1 -S2/10/0 -T15/150 -D1/1/0.5 -C0.1/0.333/0.025/0.0625 -Y1 -Zweight_test.dat -Mscak_34 -m4.5 -I1 -R-0.333/0.972/45/90/0 20090407201255351
#
# source #4 (off-diagonal):
# cap.pl -H0.02 -P1/15/60 -p1 -S2/10/0 -T15/150 -D1/1/0.5 -C0.1/0.333/0.025/0.0625 -Y1 -Zweight_test.dat -Mscak_34 -m4.5 -I1 -R0/0/90/90/90 20090407201255351
#
# source #5 (off-diagonal):
# cap.pl -H0.02 -P1/15/60 -p1 -S2/10/0 -T15/150 -D1/1/0.5 -C0.1/0.333/0.025/0.0625 -Y1 -Zweight_test.dat -Mscak_34 -m4.5 -I1 -R0/0/90/0/0 20090407201255351
#
# source #6 (off-diagonal):
# cap.pl -H0.02 -P1/15/60 -p1 -S2/10/0 -T15/150 -D1/1/0.5 -C0.1/0.333/0.025/0.0625 -Y1 -Zweight_test.dat -Mscak_34 -m4.5 -I1 -R0/0/0/90/180 20090407201255351
#
# by default, the script runs with figure generation and error checking
# turned on
import argparse
parser = argparse.ArgumentParser()
parser.add_argument('--no_checks', action='store_true')
parser.add_argument('--no_figures', action='store_true')
args = parser.parse_args()
run_checks = (not args.no_checks)
run_figures = (not args.no_figures)
from mtuq.util.cap import\
get_synthetics_cap, get_synthetics_mtuq,\
get_data_cap, compare_cap_mtuq
# the following directories correspond to the moment tensors in the list
# "sources" below
paths = []
paths += [fullpath('data/tests/benchmark_cap/20090407201255351/0')]
paths += [fullpath('data/tests/benchmark_cap/20090407201255351/1')]
paths += [fullpath('data/tests/benchmark_cap/20090407201255351/2')]
paths += [fullpath('data/tests/benchmark_cap/20090407201255351/3')]
paths += [fullpath('data/tests/benchmark_cap/20090407201255351/4')]
paths += [fullpath('data/tests/benchmark_cap/20090407201255351/5')]
paths += [fullpath('data/tests/benchmark_cap/20090407201255351/6')]
path_greens= fullpath('data/tests/benchmark_cap/greens/scak')
path_data= fullpath('data/examples/20090407201255351/*.[zrt]')
path_weights= fullpath('data/tests/benchmark_cap/20090407201255351/weights.dat')
event_id= '20090407201255351'
model= 'scak'
process_bw = ProcessData(
filter_type='Bandpass',
freq_min= 0.1,
freq_max= 0.333,
pick_type='FK_metadata',
FK_database=path_greens,
window_type='body_wave',
window_length=15.,
capuaf_file=path_weights,
)
process_sw = ProcessData(
filter_type='Bandpass',
freq_min=0.025,
freq_max=0.0625,
pick_type='FK_metadata',
FK_database=path_greens,
window_type='surface_wave',
window_length=150.,
capuaf_file=path_weights,
)
misfit_bw = Misfit(
norm='L2',
time_shift_min=-2.,
time_shift_max=0.,
time_shift_groups=['ZR'],
)
misfit_sw = Misfit(
norm='L2',
time_shift_min=-10.,
time_shift_max=0.,
time_shift_groups=['ZR','T'],
)
#
# Next we specify the source parameter grid
#
magnitude = 4.5
moment = 10.**(1.5*magnitude + 9.1) # units: N-m
sources = []
for array in [
# Mrr, Mtt, Mpp, Mrt, Mrp, Mtp
np.sqrt(1./3.)*np.array([1., 1., 1., 0., 0., 0.]), # explosion
np.array([1., 0., 0., 0., 0., 0.]), # source 1 (on-diagonal)
np.array([0., 1., 0., 0., 0., 0.]), # source 2 (on-diagonal)
np.array([0., 0., 1., 0., 0., 0.]), # source 3 (on-diagonal)
np.sqrt(1./2.)*np.array([0., 0., 0., 1., 0., 0.]), # source 4 (off-diagonal)
np.sqrt(1./2.)*np.array([0., 0., 0., 0., 1., 0.]), # source 5 (off-diagonal)
np.sqrt(1./2.)*np.array([0., 0., 0., 0., 0., 1.]), # source 6 (off-diagonal)
]:
sources += [MomentTensor(np.sqrt(2)*moment*array)]
wavelet = Trapezoid(
magnitude=magnitude)
#
# The benchmark starts now
#
print('Reading data...\n')
data = read(path_data, format='sac',
event_id=event_id,
tags=['units:m', 'type:velocity'])
data.sort_by_distance()
stations = data.get_stations()
origin = data.get_origins()[0]
print('Processing data...\n')
data_bw = data.map(process_bw)
data_sw = data.map(process_sw)
print('Reading Greens functions...\n')
db = open_db(path_greens, format='FK', model=model)
greens = db.get_greens_tensors(stations, origin)
print('Processing Greens functions...\n')
greens.convolve(wavelet)
greens_bw = greens.map(process_bw)
greens_sw = greens.map(process_sw)
depth = int(origin.depth_in_m/1000.)+1
name = '_'.join([model, str(depth), event_id])
print('Comparing waveforms...')
for _i, mt in enumerate(sources):
print(' %d of %d' % (_i+1, len(sources)))
cap_bw, cap_sw = get_synthetics_cap(
data_bw, data_sw, paths[_i], name)
mtuq_bw, mtuq_sw = get_synthetics_mtuq(
data_bw, data_sw, greens_bw, greens_sw, mt)
if run_figures:
plot_waveforms2('cap_vs_mtuq_'+str(_i)+'.png',
cap_bw, cap_sw, mtuq_bw, mtuq_sw,
stations, origin, trace_label_writer=None)
if run_checks:
compare_cap_mtuq(
cap_bw, cap_sw, mtuq_bw, mtuq_sw)
if run_figures:
# "bonus" figure comparing how CAP processes observed data with how
# MTUQ processes observed data
mtuq_sw, mtuq_bw = data_bw, data_sw
cap_sw, cap_bw = get_data_cap(
data_bw, data_sw, paths[0], name)
plot_waveforms2('cap_vs_mtuq_data.png',
cap_bw, cap_sw, mtuq_bw, mtuq_sw,
stations, origin, trace_label_writer=None, normalize=False)
print('\nSUCCESS\n')