/
wavesolver.py
373 lines (333 loc) · 15.5 KB
/
wavesolver.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
# coding: utf-8
from devito import Function, TimeFunction, warning
from devito.tools import memoized_meth
from examples.seismic.tti.operators import ForwardOperator, AdjointOperator
from examples.seismic.tti.operators import JacobianOperator, JacobianAdjOperator
from examples.seismic.tti.operators import particle_velocity_fields
from examples.checkpointing.checkpoint import DevitoCheckpoint, CheckpointOperator
from pyrevolve import Revolver
class AnisotropicWaveSolver(object):
"""
Solver object that provides operators for seismic inversion problems
and encapsulates the time and space discretization for a given problem
setup.
Parameters
----------
model : Model
Object containing the physical parameters.
geometry : AcquisitionGeometry
Geometry object that contains the source (SparseTimeFunction) and
receivers (SparseTimeFunction) and their position.
space_order : int, optional
Order of the spatial stencil discretisation. Defaults to 4.
Notes
-----
space_order must be even and it is recommended to be a multiple of 4
"""
def __init__(self, model, geometry, space_order=4, kernel='centered',
**kwargs):
self.model = model
self.model._initialize_bcs(bcs="damp")
self.geometry = geometry
self.kernel = kernel
if space_order % 2 != 0:
raise ValueError("space_order must be even but got %s"
% space_order)
if space_order % 4 != 0:
warning("It is recommended for space_order to be a multiple of 4" +
"but got %s" % space_order)
self.space_order = space_order
# Cache compiler options
self._kwargs = kwargs
@property
def dt(self):
return self.model.critical_dt
@memoized_meth
def op_fwd(self, save=False):
"""Cached operator for forward runs with buffered wavefield"""
return ForwardOperator(self.model, save=save, geometry=self.geometry,
space_order=self.space_order, kernel=self.kernel,
**self._kwargs)
@memoized_meth
def op_adj(self):
"""Cached operator for adjoint runs"""
return AdjointOperator(self.model, save=None, geometry=self.geometry,
space_order=self.space_order, kernel=self.kernel,
**self._kwargs)
@memoized_meth
def op_jac(self):
"""Cached operator for born runs"""
return JacobianOperator(self.model, save=None, geometry=self.geometry,
space_order=self.space_order, **self._kwargs)
@memoized_meth
def op_jacadj(self, save=True):
"""Cached operator for gradient runs"""
return JacobianAdjOperator(self.model, save=save, geometry=self.geometry,
space_order=self.space_order, **self._kwargs)
def forward(self, src=None, rec=None, u=None, v=None, model=None,
save=False, **kwargs):
"""
Forward modelling function that creates the necessary
data objects for running a forward modelling operator.
Parameters
----------
geometry : AcquisitionGeometry
Geometry object that contains the source (SparseTimeFunction) and
receivers (SparseTimeFunction) and their position.
u : TimeFunction, optional
The computed wavefield first component.
v : TimeFunction, optional
The computed wavefield second component.
model : Model, optional
Object containing the physical parameters.
vp : Function or float, optional
The time-constant velocity.
epsilon : Function or float, optional
The time-constant first Thomsen parameter.
delta : Function or float, optional
The time-constant second Thomsen parameter.
theta : Function or float, optional
The time-constant Dip angle (radians).
phi : Function or float, optional
The time-constant Azimuth angle (radians).
save : bool, optional
Whether or not to save the entire (unrolled) wavefield.
kernel : str, optional
Type of discretization, centered or shifted.
Returns
-------
Receiver, wavefield and performance summary.
"""
if self.kernel == 'staggered':
time_order = 1
dims = self.model.space_dimensions
stagg_u = (-dims[-1])
stagg_v = (-dims[0], -dims[1]) if self.model.grid.dim == 3 else (-dims[0])
else:
time_order = 2
stagg_u = stagg_v = None
# Source term is read-only, so re-use the default
src = src or self.geometry.src
# Create a new receiver object to store the result
rec = rec or self.geometry.rec
# Create the forward wavefield if not provided
if u is None:
u = TimeFunction(name='u', grid=self.model.grid, staggered=stagg_u,
save=self.geometry.nt if save else None,
time_order=time_order,
space_order=self.space_order)
# Create the forward wavefield if not provided
if v is None:
v = TimeFunction(name='v', grid=self.model.grid, staggered=stagg_v,
save=self.geometry.nt if save else None,
time_order=time_order,
space_order=self.space_order)
if self.kernel == 'staggered':
vx, vz, vy = particle_velocity_fields(self.model, self.space_order)
kwargs["vx"] = vx
kwargs["vz"] = vz
if vy is not None:
kwargs["vy"] = vy
model = model or self.model
# Pick vp and Thomsen parameters from model unless explicitly provided
kwargs.update(model.physical_params(**kwargs))
if self.model.dim < 3:
kwargs.pop('phi', None)
# Execute operator and return wavefield and receiver data
summary = self.op_fwd(save).apply(src=src, rec=rec, u=u, v=v,
dt=kwargs.pop('dt', self.dt), **kwargs)
return rec, u, v, summary
def adjoint(self, rec, srca=None, p=None, r=None, model=None,
save=None, **kwargs):
"""
Adjoint modelling function that creates the necessary
data objects for running an adjoint modelling operator.
Parameters
----------
geometry : AcquisitionGeometry
Geometry object that contains the source (SparseTimeFunction) and
receivers (SparseTimeFunction) and their position.
p : TimeFunction, optional
The computed wavefield first component.
r : TimeFunction, optional
The computed wavefield second component.
model : Model, optional
Object containing the physical parameters.
vp : Function or float, optional
The time-constant velocity.
epsilon : Function or float, optional
The time-constant first Thomsen parameter.
delta : Function or float, optional
The time-constant second Thomsen parameter.
theta : Function or float, optional
The time-constant Dip angle (radians).
phi : Function or float, optional
The time-constant Azimuth angle (radians).
Returns
-------
Adjoint source, wavefield and performance summary.
"""
if self.kernel == 'staggered':
time_order = 1
dims = self.model.space_dimensions
stagg_p = (-dims[-1])
stagg_r = (-dims[0], -dims[1]) if self.model.grid.dim == 3 else (-dims[0])
else:
time_order = 2
stagg_p = stagg_r = None
# Source term is read-only, so re-use the default
srca = srca or self.geometry.new_src(name='srca', src_type=None)
# Create the wavefield if not provided
if p is None:
p = TimeFunction(name='p', grid=self.model.grid, staggered=stagg_p,
time_order=time_order,
space_order=self.space_order)
# Create the wavefield if not provided
if r is None:
r = TimeFunction(name='r', grid=self.model.grid, staggered=stagg_r,
time_order=time_order,
space_order=self.space_order)
if self.kernel == 'staggered':
vx, vz, vy = particle_velocity_fields(self.model, self.space_order)
kwargs["vx"] = vx
kwargs["vz"] = vz
if vy is not None:
kwargs["vy"] = vy
model = model or self.model
# Pick vp and Thomsen parameters from model unless explicitly provided
kwargs.update(model.physical_params(**kwargs))
if self.model.dim < 3:
kwargs.pop('phi', None)
# Execute operator and return wavefield and receiver data
summary = self.op_adj().apply(srca=srca, rec=rec, p=p, r=r,
dt=kwargs.pop('dt', self.dt),
time_m=0 if time_order == 1 else None,
**kwargs)
return srca, p, r, summary
def jacobian(self, dm, src=None, rec=None, u0=None, v0=None, du=None, dv=None,
model=None, save=None, kernel='centered', **kwargs):
"""
Linearized Born modelling function that creates the necessary
data objects for running an adjoint modelling operator.
Parameters
----------
src : SparseTimeFunction or array_like, optional
Time series data for the injected source term.
rec : SparseTimeFunction or array_like, optional
The interpolated receiver data.
u : TimeFunction, optional
The computed background wavefield first component.
v : TimeFunction, optional
The computed background wavefield second component.
du : TimeFunction, optional
The computed perturbed wavefield first component.
dv : TimeFunction, optional
The computed perturbed wavefield second component.
model : Model, optional
Object containing the physical parameters.
vp : Function or float, optional
The time-constant velocity.
epsilon : Function or float, optional
The time-constant first Thomsen parameter.
delta : Function or float, optional
The time-constant second Thomsen parameter.
theta : Function or float, optional
The time-constant Dip angle (radians).
phi : Function or float, optional
The time-constant Azimuth angle (radians).
"""
if kernel != 'centered':
raise ValueError('Only centered kernel is supported for the jacobian')
dt = kwargs.pop('dt', self.dt)
# Source term is read-only, so re-use the default
src = src or self.geometry.src
# Create a new receiver object to store the result
rec = rec or self.geometry.rec
# Create the forward wavefields u, v du and dv if not provided
u0 = u0 or TimeFunction(name='u0', grid=self.model.grid,
time_order=2, space_order=self.space_order)
v0 = v0 or TimeFunction(name='v0', grid=self.model.grid,
time_order=2, space_order=self.space_order)
du = du or TimeFunction(name='du', grid=self.model.grid,
time_order=2, space_order=self.space_order)
dv = dv or TimeFunction(name='dv', grid=self.model.grid,
time_order=2, space_order=self.space_order)
model = model or self.model
# Pick vp and Thomsen parameters from model unless explicitly provided
kwargs.update(model.physical_params(**kwargs))
if self.model.dim < 3:
kwargs.pop('phi', None)
# Execute operator and return wavefield and receiver data
summary = self.op_jac().apply(dm=dm, u0=u0, v0=v0, du=du, dv=dv, src=src,
rec=rec, dt=dt, **kwargs)
return rec, u0, v0, du, dv, summary
def jacobian_adjoint(self, rec, u0, v0, du=None, dv=None, dm=None, model=None,
checkpointing=False, kernel='centered', **kwargs):
"""
Gradient modelling function for computing the adjoint of the
Linearized Born modelling function, ie. the action of the
Jacobian adjoint on an input data.
Parameters
----------
rec : SparseTimeFunction
Receiver data.
u0 : TimeFunction
The computed background wavefield.
v0 : TimeFunction, optional
The computed background wavefield.
du : Function or float
The computed perturbed wavefield.
dv : Function or float
The computed perturbed wavefield.
dm : Function, optional
Stores the gradient field.
model : Model, optional
Object containing the physical parameters.
vp : Function or float, optional
The time-constant velocity.
epsilon : Function or float, optional
The time-constant first Thomsen parameter.
delta : Function or float, optional
The time-constant second Thomsen parameter.
theta : Function or float, optional
The time-constant Dip angle (radians).
phi : Function or float, optional
The time-constant Azimuth angle (radians).
Returns
-------
Gradient field and performance summary.
"""
if kernel != 'centered':
raise ValueError('Only centered kernel is supported for the jacobian_adj')
dt = kwargs.pop('dt', self.dt)
# Gradient symbol
dm = dm or Function(name='dm', grid=self.model.grid)
# Create the perturbation wavefields if not provided
du = du or TimeFunction(name='du', grid=self.model.grid,
time_order=2, space_order=self.space_order)
dv = dv or TimeFunction(name='dv', grid=self.model.grid,
time_order=2, space_order=self.space_order)
model = model or self.model
# Pick vp and Thomsen parameters from model unless explicitly provided
kwargs.update(model.physical_params(**kwargs))
if self.model.dim < 3:
kwargs.pop('phi', None)
if checkpointing:
u0 = TimeFunction(name='u0', grid=self.model.grid,
time_order=2, space_order=self.space_order)
v0 = TimeFunction(name='v0', grid=self.model.grid,
time_order=2, space_order=self.space_order)
cp = DevitoCheckpoint([u0, v0])
n_checkpoints = None
wrap_fw = CheckpointOperator(self.op_fwd(save=False), src=self.geometry.src,
u=u0, v=v0, dt=dt, **kwargs)
wrap_rev = CheckpointOperator(self.op_jacadj(save=False), u0=u0, v0=v0,
du=du, dv=dv, rec=rec, dm=dm, dt=dt, **kwargs)
# Run forward
wrp = Revolver(cp, wrap_fw, wrap_rev, n_checkpoints, rec.data.shape[0]-2)
wrp.apply_forward()
summary = wrp.apply_reverse()
else:
summary = self.op_jacadj().apply(rec=rec, dm=dm, u0=u0, v0=v0, du=du, dv=dv,
dt=dt, **kwargs)
return dm, summary