forked from fpicetti/occamypy
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpreset_models.py
358 lines (299 loc) · 15.4 KB
/
preset_models.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
import os
import numpy as np
from .model import SeismicModel
__all__ = ['demo_model']
def demo_model(preset, **kwargs):
"""
Utility function to create preset `Model` objects for
demonstration and testing purposes. The particular presets are ::
* `constant-isotropic` : Constant velocity (1.5 km/sec) isotropic model
* `constant-tti` : Constant anisotropic model. Velocity is 1.5 km/sec and
Thomsen parameters are epsilon=.3, delta=.2, theta = .7rad
and phi=.35rad for 3D. 2d/3d is defined from the input shape
* 'layers-isotropic': Simple n-layered model with velocities ranging from 1.5 km/s
to 3.5 km/s in the top and bottom layer respectively.
2d/3d is defined from the input shape
* 'layers-elastic': Simple n-layered model with velocities ranging from 1.5 km/s
to 3.5 km/s in the top and bottom layer respectively.
Vs is set to .5 vp and 0 in the top layer.
* 'layers-viscoelastic': Simple two layers viscoelastic model.
* 'layers-tti': Simple n-layered model with velocities ranging from 1.5 km/s
to 3.5 km/s in the top and bottom layer respectively.
Thomsen parameters in the top layer are 0 and in the lower layers
are scaled versions of vp.
2d/3d is defined from the input shape
* 'circle-isotropic': Simple camembert model with velocities 1.5 km/s
and 2.5 km/s in a circle at the center. 2D only.
* 'marmousi2d-isotropic': Loads the 2D Marmousi data set from the given
filepath. Requires the ``devitocodes/data`` repository
to be available on your machine.
* 'marmousi2d-tti': Loads the 2D Marmousi data set from the given
filepath. Requires the ``devitocodes/data`` repository
to be available on your machine.
* 'marmousi3d-tti': Loads the 2D Marmousi data set from the given
filepath. Requires the ``devitocodes/data`` repository
to be available on your machine.
"""
space_order = kwargs.pop('space_order', 2)
shape = kwargs.pop('shape', (101, 101))
spacing = kwargs.pop('spacing', tuple([10. for _ in shape]))
origin = kwargs.pop('origin', tuple([0. for _ in shape]))
nbl = kwargs.pop('nbl', 10)
dtype = kwargs.pop('dtype', np.float32)
vp = kwargs.pop('vp', 1.5)
nlayers = kwargs.pop('nlayers', 3)
fs = kwargs.pop('fs', False)
if preset.lower() in ['constant-elastic']:
# op constant single-layer model in a 2D or 3D domain
# with velocity 1.5 km/s.
vs = 0.5 * vp
b = 1.0
return SeismicModel(space_order=space_order, vp=vp, vs=vs, b=b,
origin=origin, shape=shape, dtype=dtype, spacing=spacing,
nbl=nbl, **kwargs)
if preset.lower() in ['constant-viscoelastic']:
# op constant single-layer model in a 2D or 3D domain
# with velocity 2.2 km/s.
qp = kwargs.pop('qp', 100.)
vs = kwargs.pop('vs', 1.2)
qs = kwargs.pop('qs', 70.)
b = 1/2.
return SeismicModel(space_order=space_order, vp=vp, qp=qp, vs=vs,
qs=qs, b=b, origin=origin, shape=shape,
dtype=dtype, spacing=spacing, nbl=nbl,
**kwargs)
if preset.lower() in ['constant-isotropic']:
# op constant single-layer model in a 2D or 3D domain
# with velocity 1.5 km/s.
return SeismicModel(space_order=space_order, vp=vp, origin=origin, shape=shape,
dtype=dtype, spacing=spacing, nbl=nbl, fs=fs, **kwargs)
if preset.lower() in ['constant-viscoacoustic']:
# op constant single-layer model in a 2D or 3D domain
# with velocity 1.5 km/s.
qp = kwargs.pop('qp', 100.)
b = 1/2.
return SeismicModel(space_order=space_order, vp=vp, qp=qp, b=b, nbl=nbl,
dtype=dtype, origin=origin, shape=shape,
spacing=spacing, **kwargs)
elif preset.lower() in ['constant-tti']:
# op constant single-layer model in a 2D or 3D domain
# with velocity 1.5 km/s.
v = np.empty(shape, dtype=dtype)
v[:] = 1.5
epsilon = .3*np.ones(shape, dtype=dtype)
delta = .2*np.ones(shape, dtype=dtype)
theta = .7*np.ones(shape, dtype=dtype)
phi = None
if len(shape) > 2:
phi = .35*np.ones(shape, dtype=dtype)
return SeismicModel(space_order=space_order, vp=v, origin=origin, shape=shape,
dtype=dtype, spacing=spacing, nbl=nbl, epsilon=epsilon,
delta=delta, theta=theta, phi=phi, bcs="damp", **kwargs)
elif preset.lower() in ['layers-isotropic']:
# op n-layers model in a 2D or 3D domain with two different
# velocities split across the height dimension:
# By default, the top part of the domain has 1.5 km/s,
# and the bottom part of the domain has 2.5 km/s.
vp_top = kwargs.pop('vp_top', 1.5)
vp_bottom = kwargs.pop('vp_bottom', 3.5)
# Define a velocity profile in km/s
v = np.empty(shape, dtype=dtype)
v[:] = vp_top # Top velocity (background)
vp_i = np.linspace(vp_top, vp_bottom, nlayers)
for i in range(1, nlayers):
v[..., i*int(shape[-1] / nlayers):] = vp_i[i] # Bottom velocity
return SeismicModel(space_order=space_order, vp=v, origin=origin, shape=shape,
dtype=dtype, spacing=spacing, nbl=nbl, bcs="damp",
fs=fs, **kwargs)
elif preset.lower() in ['layers-elastic']:
# op n-layers model in a 2D or 3D domain with two different
# velocities split across the height dimension:
# By default, the top part of the domain has 1.5 km/s,
# and the bottom part of the domain has 2.5 km/s.
vp_top = kwargs.pop('vp_top', 1.5)
vp_bottom = kwargs.pop('vp_bottom', 3.5)
# Define a velocity profile in km/s
v = np.empty(shape, dtype=dtype)
v[:] = vp_top # Top velocity (background)
vp_i = np.linspace(vp_top, vp_bottom, nlayers)
for i in range(1, nlayers):
v[..., i*int(shape[-1] / nlayers):] = vp_i[i] # Bottom velocity
vs = 0.5 * v[:]
b = 1 / (0.31 * (1e3*v)**0.25)
b[v < 1.51] = 1.0
vs[v < 1.51] = 0.0
return SeismicModel(space_order=space_order, vp=v, vs=vs, b=b,
origin=origin, shape=shape,
dtype=dtype, spacing=spacing, nbl=nbl, **kwargs)
elif preset.lower() in ['layers-viscoelastic', 'twolayer-viscoelastic',
'2layer-viscoelastic']:
# op two-layer model in a 2D or 3D domain with two different
# velocities split across the height dimension:
# By default, the top part of the domain has 1.6 km/s,
# and the bottom part of the domain has 2.2 km/s.
ratio = kwargs.pop('ratio', 3)
vp_top = kwargs.pop('vp_top', 1.6)
qp_top = kwargs.pop('qp_top', 40.)
vs_top = kwargs.pop('vs_top', 0.4)
qs_top = kwargs.pop('qs_top', 30.)
b_top = kwargs.pop('b_top', 1/1.3)
vp_bottom = kwargs.pop('vp_bottom', 2.2)
qp_bottom = kwargs.pop('qp_bottom', 100.)
vs_bottom = kwargs.pop('vs_bottom', 1.2)
qs_bottom = kwargs.pop('qs_bottom', 70.)
b_bottom = kwargs.pop('b_bottom', 1/2.)
# Define a velocity profile in km/s
vp = np.empty(shape, dtype=dtype)
qp = np.empty(shape, dtype=dtype)
vs = np.empty(shape, dtype=dtype)
qs = np.empty(shape, dtype=dtype)
b = np.empty(shape, dtype=dtype)
# Top and bottom P-wave velocity
vp[:] = vp_top
vp[..., int(shape[-1] / ratio):] = vp_bottom
# Top and bottom P-wave quality factor
qp[:] = qp_top
qp[..., int(shape[-1] / ratio):] = qp_bottom
# Top and bottom S-wave velocity
vs[:] = vs_top
vs[..., int(shape[-1] / ratio):] = vs_bottom
# Top and bottom S-wave quality factor
qs[:] = qs_top
qs[..., int(shape[-1] / ratio):] = qs_bottom
# Top and bottom density
b[:] = b_top
b[..., int(shape[-1] / ratio):] = b_bottom
return SeismicModel(space_order=space_order, vp=vp, qp=qp,
vs=vs, qs=qs, b=b, origin=origin,
shape=shape, dtype=dtype, spacing=spacing,
nbl=nbl, **kwargs)
elif preset.lower() in ['layers-tti', 'layers-tti-noazimuth']:
# op n-layers model in a 2D or 3D domain with two different
# velocities split across the height dimension:
# By default, the top part of the domain has 1.5 km/s,
# and the bottom part of the domain has 2.5 km/s.\
vp_top = kwargs.pop('vp_top', 1.5)
vp_bottom = kwargs.pop('vp_bottom', 3.5)
# Define a velocity profile in km/s
v = np.empty(shape, dtype=dtype)
v[:] = vp_top # Top velocity (background)
vp_i = np.linspace(vp_top, vp_bottom, nlayers)
for i in range(1, nlayers):
v[..., i*int(shape[-1] / nlayers):] = vp_i[i] # Bottom velocity
epsilon = .3*(v - 1.5)
delta = .2*(v - 1.5)
theta = .5*(v - 1.5)
phi = None
if len(shape) > 2 and preset.lower() not in ['layers-tti-noazimuth']:
phi = .25*(v - 1.5)
model = SeismicModel(space_order=space_order, vp=v, origin=origin, shape=shape,
dtype=dtype, spacing=spacing, nbl=nbl, epsilon=epsilon,
delta=delta, theta=theta, phi=phi, bcs="damp", **kwargs)
if kwargs.get('smooth', False):
if len(shape) > 2 and preset.lower() not in ['layers-tti-noazimuth']:
model.smooth(('epsilon', 'delta', 'theta', 'phi'))
else:
model.smooth(('epsilon', 'delta', 'theta'))
return model
elif preset.lower() in ['circle-isotropic']:
# op simple circle in a 2D domain with a background velocity.
# By default, the circle velocity is 2.5 km/s,
# and the background veloity is 3.0 km/s.
vp = kwargs.pop('vp_circle', 3.0)
vp_background = kwargs.pop('vp_background', 2.5)
r = kwargs.pop('r', 15)
# Only a 2D preset is available currently
assert(len(shape) == 2)
v = np.empty(shape, dtype=dtype)
v[:] = vp_background
a, b = shape[0] / 2, shape[1] / 2
y, x = np.ogrid[-a:shape[0]-a, -b:shape[1]-b]
v[x*x + y*y <= r*r] = vp
return SeismicModel(space_order=space_order, vp=v, origin=origin, shape=shape,
dtype=dtype, spacing=spacing, nbl=nbl, bcs="damp",
fs=fs, **kwargs)
elif preset.lower() in ['marmousi-isotropic', 'marmousi2d-isotropic']:
shape = (1601, 401)
spacing = (7.5, 7.5)
origin = (0., 0.)
nbl = kwargs.pop('nbl', 20)
# Read 2D Marmousi model from devitocodes/data repo
data_path = kwargs.get('data_path', None)
if data_path is None:
raise ValueError("Path to devitocodes/data not found! Please specify with "
"'data_path=<path/to/devitocodes/data>'")
path = os.path.join(data_path, 'Simple2D/vp_marmousi_bi')
v = np.fromfile(path, dtype='float32', sep="")
v = v.reshape(shape)
# Cut the model to make it slightly cheaper
v = v[301:-300, :]
return SeismicModel(space_order=space_order, vp=v, origin=origin, shape=v.shape,
dtype=np.float32, spacing=spacing, nbl=nbl, bcs="damp",
fs=fs, **kwargs)
elif preset.lower() in ['marmousi-tti2d', 'marmousi2d-tti',
'marmousi-tti3d', 'marmousi3d-tti']:
shape = (201, 201, 70)
nbl = kwargs.pop('nbl', 20)
# Read 2D Marmousi model from devitocodes/data repo
data_path = kwargs.pop('data_path', None)
if data_path is None:
raise ValueError("Path to devitocodes/data not found! Please specify with "
"'data_path=<path/to/devitocodes/data>'")
path = os.path.join(data_path, 'marmousi3D/vp_marmousi_bi')
# velocity
vp = 1e-3 * np.fromfile(os.path.join(data_path, 'marmousi3D/MarmousiVP.raw'),
dtype='float32', sep="")
vp = vp.reshape(shape)
# Epsilon, in % in file, resale between 0 and 1
epsilon = np.fromfile(os.path.join(data_path, 'marmousi3D/MarmousiEps.raw'),
dtype='float32', sep="") * 1e-2
epsilon = epsilon.reshape(shape)
# Delta, in % in file, resale between 0 and 1
delta = np.fromfile(os.path.join(data_path, 'marmousi3D/MarmousiDelta.raw'),
dtype='float32', sep="") * 1e-2
delta = delta.reshape(shape)
# Theta, in degrees in file, resale in radian
theta = np.fromfile(os.path.join(data_path, 'marmousi3D/MarmousiTilt.raw'),
dtype='float32', sep="")
theta = np.float32(np.pi / 180 * theta.reshape(shape))
if preset.lower() in ['marmousi-tti3d', 'marmousi3d-tti']:
# Phi, in degrees in file, resale in radian
phi = np.fromfile(os.path.join(data_path, 'marmousi3D/Azimuth.raw'),
dtype='float32', sep="")
phi = np.float32(np.pi / 180 * phi.reshape(shape))
else:
vp = vp[101, :, :]
epsilon = epsilon[101, :, :]
delta = delta[101, :, :]
theta = theta[101, :, :]
shape = vp.shape
phi = None
spacing = tuple([10.0]*len(shape))
origin = tuple([0.0]*len(shape))
return SeismicModel(space_order=space_order, vp=vp, origin=origin, shape=shape,
dtype=np.float32, spacing=spacing, nbl=nbl, epsilon=epsilon,
delta=delta, theta=theta, phi=phi, bcs="damp", **kwargs)
elif preset.lower() in ['layers-viscoacoustic']:
# op n-layers model in a 2D or 3D domain with two different
# velocities split across the height dimension:
# By default, the top part of the domain has 1.5 km/s,
# and the bottom part of the domain has 3.5 km/s.
# Define a velocity profile in km/s
vp = np.empty(shape, dtype=dtype)
qp = np.empty(shape, dtype=dtype)
# Top and bottom P-wave velocity
vp_top = kwargs.pop('vp_top', 1.5)
vp_bottom = kwargs.pop('vp_bottom', 3.5)
# Define a velocity profile in km/s
vp = np.empty(shape, dtype=dtype)
vp[:] = vp_top # Top velocity (background)
vp_i = np.linspace(vp_top, vp_bottom, nlayers)
for i in range(1, nlayers):
vp[..., i*int(shape[-1] / nlayers):] = vp_i[i] # Bottom velocity
qp[:] = 3.516*((vp[:]*1000.)**2.2)*10**(-6) # Li's empirical formula
b = 1 / (0.31*(vp[:]*1000.)**0.25) # Gardner's relation
return SeismicModel(space_order=space_order, vp=vp, qp=qp, b=b, nbl=nbl,
dtype=dtype, origin=origin, shape=shape,
spacing=spacing, **kwargs)
else:
raise ValueError("Unknown model preset name")