-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathGridSearch.FullMomentTensor_SA.py
222 lines (159 loc) · 5.48 KB
/
GridSearch.FullMomentTensor_SA.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
#!/usr/bin/env python
import os
import numpy as np
from mtuq import read, open_db, download_greens_tensors
from mtuq.event import Origin
from mtuq.graphics import plot_data_greens2, plot_beachball, plot_misfit_lune
from mtuq.grid import FullMomentTensorGridSemiregular
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__':
#
# Carries out grid search over all moment tensor parameters
#
# USAGE
# mpirun -n <NPROC> python GridSearch.FullMomentTensor.py
#
path_data= fullpath('/home/jovyan/scoped/pysep/20140828081339000/*.[zrt]')
path_weights= fullpath('/home/jovyan/scoped/pysep/20140828081339000/weights.dat')
event_id= '20140828081339000'
model= 'ak135'
#
# Body and surface wave measurements will be made separately
#
process_bw = ProcessData(
filter_type='Bandpass',
freq_min= 0.1,
freq_max= 0.333,
pick_type='taup',
taup_model=model,
window_type='body_wave',
window_length=15.,
capuaf_file=path_weights,
)
process_sw = ProcessData(
filter_type='Bandpass',
freq_min=0.02,
freq_max=0.05,
pick_type='taup',
taup_model=model,
window_type='surface_wave',
window_length=150.,
capuaf_file=path_weights,
)
#
# For our objective function, we will use a sum of body and surface wave
# contributions
#
misfit_bw = Misfit(
norm='L2',
time_shift_min=-4.,
time_shift_max=+4.,
time_shift_groups=['ZR'],
)
misfit_sw = Misfit(
norm='L2',
time_shift_min=-15.,
time_shift_max=+15.,
time_shift_groups=['ZR','T'],
)
#
# User-supplied weights control how much each station contributes to the
# objective function
#
station_id_list = parse_station_codes(path_weights)
#
# Next, we specify the moment tensor grid and source-time function
#
grid = FullMomentTensorGridSemiregular(
npts_per_axis=15,
magnitudes=[5.3])
wavelet = Trapezoid(
magnitude=5.3)
#
# Origin time and location will be fixed. For an example in which they
# vary, see examples/GridSearch.DoubleCouple+Magnitude+Depth.py
#
# See also Dataset.get_origins(), which attempts to create Origin objects
# from waveform metadata
#
origin = Origin({
'time': '2014-08-28T08:13:39.00000Z',
'latitude': 64.654,
'longitude': -17.385,
'depth_in_m': 7000.0,
})
from mpi4py import MPI
comm = MPI.COMM_WORLD
#
# The main I/O work starts now
#
if comm.rank==0:
print('Reading data...\n')
data = read(path_data, format='sac',
event_id=event_id,
station_id_list=station_id_list,
tags=['units:cm', 'type:velocity'])
data.sort_by_distance()
stations = data.get_stations()
print('Processing data...\n')
data_bw = data.map(process_bw)
data_sw = data.map(process_sw)
print('Reading Greens functions...\n')
greens = download_greens_tensors(stations, origin, model)
print('Processing Greens functions...\n')
greens.convolve(wavelet)
greens_bw = greens.map(process_bw)
greens_sw = greens.map(process_sw)
else:
stations = None
data_bw = None
data_sw = None
greens_bw = None
greens_sw = None
stations = comm.bcast(stations, root=0)
data_bw = comm.bcast(data_bw, root=0)
data_sw = comm.bcast(data_sw, root=0)
greens_bw = comm.bcast(greens_bw, root=0)
greens_sw = comm.bcast(greens_sw, root=0)
#
# The main computational work starts now
#
if comm.rank==0:
print('Evaluating body wave misfit...\n')
results_bw = grid_search(
data_bw, greens_bw, misfit_bw, origin, grid)
if comm.rank==0:
print('Evaluating surface wave misfit...\n')
results_sw = grid_search(
data_sw, greens_sw, misfit_sw, origin, grid)
if comm.rank==0:
results = results_bw + results_sw
# array index corresponding to minimum misfit
idx = results.idxmin('source')
best_source = grid.get(idx)
lune_dict = grid.get_dict(idx)
mt_dict = grid.get(idx).as_dict()
#
# Generate figures and save results
#
print('Generating figures...\n')
plot_data_greens2(event_id+'FMT_waveforms.png',
data_bw, data_sw, greens_bw, greens_sw, process_bw, process_sw,
misfit_bw, misfit_sw, stations, origin, best_source, lune_dict)
plot_beachball(event_id+'FMT_beachball.png',
best_source, stations, origin)
plot_misfit_lune(event_id+'FMT_misfit.png', results)
plot_misfit_lune(event_id+'FMT_misfit_mt.png', results, show_mt=True)
plot_misfit_lune(event_id+'FMT_misfit_tradeoff.png', results, show_tradeoffs=True)
print('Saving results...\n')
merged_dict = merge_dicts(lune_dict, mt_dict, origin,
{'M0': best_source.moment(), 'Mw': best_source.magnitude()})
# save best-fitting source
save_json(event_id+'FMT_solution.json', merged_dict)
# save misfit surface
results.save(event_id+'FMT_misfit.nc')
print('\nFinished\n')