-
Notifications
You must be signed in to change notification settings - Fork 220
/
operators.py
executable file
·652 lines (491 loc) · 19 KB
/
operators.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
import sympy as sp
import numpy as np
from devito import (Eq, Operator, VectorTimeFunction, TimeFunction, Function, NODE,
div, grad, Inc)
from examples.seismic import PointSource, Receiver
def src_rec(p, model, geometry, **kwargs):
"""
Forward case: Source injection and receiver interpolation
Adjoint case: Receiver injection and source interpolation
"""
dt = model.grid.time_dim.spacing
m = model.m
# Source symbol with input wavelet
src = PointSource(name="src", grid=model.grid, time_range=geometry.time_axis,
npoint=geometry.nsrc)
rec = Receiver(name='rec', grid=model.grid, time_range=geometry.time_axis,
npoint=geometry.nrec)
forward = kwargs.get('forward', True)
time_order = p.time_order
if forward:
# The source injection term
if(time_order == 1):
src_term = src.inject(field=p.forward, expr=src * dt)
else:
src_term = src.inject(field=p.forward, expr=src * dt**2 / m)
# Create interpolation expression for receivers
rec_term = rec.interpolate(expr=p)
else:
# Construct expression to inject receiver values
if(time_order == 1):
rec_term = rec.inject(field=p.backward, expr=rec * dt)
else:
rec_term = rec.inject(field=p.backward, expr=rec * dt**2 / m)
# Create interpolation expression for the adjoint-source
src_term = src.interpolate(expr=p)
return src_term + rec_term
def sls_1st_order(model, geometry, p, **kwargs):
"""
Implementation of the 1st order viscoacoustic wave-equation
from Blanch and Symes (1995) / Dutta and Schuster (2014).
https://library.seg.org/doi/pdf/10.1190/1.1822695
https://library.seg.org/doi/pdf/10.1190/geo2013-0414.1
Parameters
----------
p : TimeFunction
Pressure field.
"""
forward = kwargs.get('forward', True)
space_order = p.space_order
save = kwargs.get('save', False)
save_t = geometry.nt if save else None
s = model.grid.stepping_dim.spacing
b = model.b
vp = model.vp
damp = model.damp
qp = model.qp
f0 = geometry._f0
# Particle Velocity
v = kwargs.pop('v')
# The stress relaxation parameter
t_s = (sp.sqrt(1.+1./qp**2)-1./qp)/f0
# The strain relaxation parameter
t_ep = 1./(f0**2*t_s)
# The relaxation time
tt = (t_ep/t_s)-1.
# Density
rho = 1. / b
# Bulk modulus
bm = rho * vp**2
# Attenuation Memory variable.
r = TimeFunction(name="r", grid=model.grid, time_order=1, space_order=space_order,
save=save_t, staggered=NODE)
if forward:
# Define PDE
pde_v = v - s * b * grad(p)
u_v = Eq(v.forward, damp * pde_v)
pde_r = r - s * (1. / t_s) * r - s * (1. / t_s) * tt * bm * div(v.forward)
u_r = Eq(r.forward, damp * pde_r)
pde_p = p - s * bm * (tt + 1.) * div(v.forward) - s * r.forward
u_p = Eq(p.forward, damp * pde_p)
return [u_v, u_r, u_p]
else:
# Define PDE
pde_r = r - s * (1. / t_s) * r - s * p
u_r = Eq(r.backward, damp * pde_r)
pde_v = v + s * grad(bm * (1. + tt) * p) + s * \
grad((1. / t_s) * bm * tt * r.backward)
u_v = Eq(v.backward, damp * pde_v)
pde_p = p + s * div(b * v.backward)
u_p = Eq(p.backward, damp * pde_p)
return [u_r, u_v, u_p]
def sls_2nd_order(model, geometry, p, r=None, **kwargs):
"""
Implementation of the 2nd order viscoacoustic wave-equation from Bai (2014).
https://library.seg.org/doi/10.1190/geo2013-0030.1
Parameters
----------
p : TimeFunction
Pressure field.
r : TimeFunction
Attenuation Memory variable.
"""
forward = kwargs.get('forward', True)
space_order = p.space_order
save = kwargs.get('save', False)
save_t = geometry.nt if save else None
s = model.grid.stepping_dim.spacing
b = model.b
vp = model.vp
damp = model.damp
qp = model.qp
f0 = geometry._f0
q = kwargs.get('q', 0)
# The stress relaxation parameter
t_s = (sp.sqrt(1.+1./qp**2)-1./qp)/f0
# The strain relaxation parameter
t_ep = 1./(f0**2*t_s)
# The relaxation time
tt = (t_ep/t_s)-1.
# Density
rho = 1. / b
# Bulk modulus
bm = rho * vp**2
# Attenuation Memory variable.
r = r or TimeFunction(name="r", grid=model.grid, time_order=2,
space_order=space_order, save=save_t, staggered=NODE)
if forward:
pde_r = r + s * (tt / t_s) * rho * div(b * grad(p, shift=.5), shift=-.5) - \
s * (1. / t_s) * r
u_r = Eq(r.forward, damp * pde_r)
pde_p = 2. * p - damp * p.backward + s**2 * bm * (1. + tt) * \
div(b * grad(p, shift=.5), shift=-.5) - s**2 * vp**2 * \
r.forward + s**2 * vp**2 * q
u_p = Eq(p.forward, damp * pde_p)
return [u_r, u_p]
else:
pde_r = r + s * (tt / t_s) * p - s * (1. / t_s) * r
u_r = Eq(r.backward, damp * pde_r)
pde_p = 2. * p - damp * p.forward + s**2 * vp**2 * \
div(b * grad((1. + tt) * rho * p, shift=.5), shift=-.5) - s**2 * vp**2 * \
div(b * grad(rho * r.backward, shift=.5), shift=-.5)
u_p = Eq(p.backward, damp * pde_p)
return [u_r, u_p]
def ren_1st_order(model, geometry, p, **kwargs):
"""
Implementation of the 1st order viscoacoustic wave-equation from Ren et al. (2014).
https://academic.oup.com/gji/article/197/2/948/616510
Parameters
----------
p : TimeFunction
Pressure field.
"""
forward = kwargs.get('forward', True)
s = model.grid.stepping_dim.spacing
f0 = geometry._f0
vp = model.vp
b = model.b
qp = model.qp
damp = model.damp
# Particle velocity
v = kwargs.pop('v')
# Angular frequency
w0 = 2. * np.pi * f0
# Density
rho = 1. / b
eta = vp**2 / (w0 * qp)
# Bulk modulus
bm = rho * vp**2
if forward:
# Define PDE
pde_v = v - s * b * grad(p)
u_v = Eq(v.forward, damp * pde_v)
pde_p = p - s * bm * div(v.forward) + \
s * ((vp**2 * rho) / (w0 * qp)) * div(b * grad(p, shift=.5), shift=-.5)
u_p = Eq(p.forward, damp * pde_p)
return [u_v, u_p]
else:
pde_v = v + s * grad(bm * p)
u_v = Eq(v.backward, pde_v * damp)
pde_p = p + s * div(b * grad(rho * eta * p, shift=.5), shift=-.5) + \
s * div(b * v.backward)
u_p = Eq(p.backward, pde_p * damp)
return [u_v, u_p]
def ren_2nd_order(model, geometry, p, **kwargs):
"""
Implementation of the 2nd order viscoacoustic wave-equation from Ren et al. (2014).
https://library.seg.org/doi/pdf/10.1190/1.2714334
Parameters
----------
p : TimeFunction
Pressure field.
"""
forward = kwargs.get('forward', True)
s = model.grid.stepping_dim.spacing
f0 = geometry._f0
vp = model.vp
b = model.b
qp = model.qp
damp = model.damp
# Angular frequency
w0 = 2. * np.pi * f0
# Density
rho = 1. / b
eta = vp**2 / (w0 * qp)
# Bulk modulus
bm = rho * vp**2
if forward:
pde_p = 2. * p - damp * p.backward + s**2 * bm * \
div(b * grad(p, shift=.5), shift=-.5) + s**2 * eta * rho * \
div(b * grad(p - p.backward, shift=.5) / s, shift=-.5)
u_p = Eq(p.forward, damp * pde_p)
return [u_p]
else:
pde_p = 2. * p - damp * p.forward + s**2 * \
div(b * grad(bm * p, shift=.5), shift=-.5) - s**2 * \
div(b * grad(((p.forward - p) / s) * rho * eta, shift=.5), shift=-.5)
u_p = Eq(p.backward, damp * pde_p)
return [u_p]
def deng_1st_order(model, geometry, p, **kwargs):
"""
Implementation of the 1st order viscoacoustic wave-equation
from Deng and McMechan (2007).
https://library.seg.org/doi/pdf/10.1190/1.2714334
Parameters
----------
p : TimeFunction
Pressure field.
"""
forward = kwargs.get('forward', True)
s = model.grid.stepping_dim.spacing
f0 = geometry._f0
vp = model.vp
b = model.b
qp = model.qp
damp = model.damp
# Particle velocity
v = kwargs.pop('v')
# Angular frequency
w0 = 2. * np.pi * f0
# Density
rho = 1. / b
# Bulk modulus
bm = rho * vp**2
if forward:
# Define PDE
pde_v = v - s * b * grad(p)
u_v = Eq(v.forward, damp * pde_v)
pde_p = p - s * bm * div(v.forward) - s * (w0 / qp) * p
u_p = Eq(p.forward, damp * pde_p)
return [u_v, u_p]
else:
pde_v = v + s * grad(bm * p)
u_v = Eq(v.backward, pde_v * damp)
pde_p = p + s * div(b * v.backward) - s * (w0 / qp) * p
u_p = Eq(p.backward, pde_p * damp)
return [u_v, u_p]
def deng_2nd_order(model, geometry, p, **kwargs):
"""
Implementation of the 2nd order viscoacoustic wave-equation
from Deng and McMechan (2007).
https://library.seg.org/doi/pdf/10.1190/1.2714334
Parameters
----------
p : TimeFunction
Pressure field.
"""
forward = kwargs.get('forward', True)
s = model.grid.stepping_dim.spacing
f0 = geometry._f0
vp = model.vp
b = model.b
qp = model.qp
damp = model.damp
# Angular frequency
w0 = 2. * np.pi * f0
# Density
rho = 1. / b
bm = rho * vp**2
if forward:
pde_p = 2. * p - damp*p.backward + s**2 * bm * \
div(b * grad(p, shift=.5), shift=-.5) - s**2 * w0/qp * (p - p.backward)/s
u_p = Eq(p.forward, damp * pde_p)
return [u_p]
else:
pde_p = 2. * p - damp * p.forward + s**2 * w0 / qp * (p.forward - p) / s + \
s * s * div(b * grad(bm * p, shift=.5), shift=-.5)
u_p = Eq(p.backward, damp * pde_p)
return [u_p]
def sls(model, geometry, p, forward=True, **kwargs):
"""
Implementation of the 1st order viscoacoustic wave-equation
from Blanch and Symes (1995) / Dutta and Schuster (2014) and
Implementation of the 2nd order viscoacoustic wave-equation from Bai (2014).
https://library.seg.org/doi/pdf/10.1190/1.1822695
https://library.seg.org/doi/pdf/10.1190/geo2013-0414.1
https://library.seg.org/doi/10.1190/geo2013-0030.1
Parameters
----------
p : TimeFunction
Pressure field.
"""
time_order = p.time_order
eq_stencil = stencils[('sls', time_order)]
eqn = eq_stencil(model, geometry, p, forward=forward, **kwargs)
return eqn
def ren(model, geometry, p, forward=True, **kwargs):
"""
Implementation of the 1st and 2nd order viscoacoustic wave-equation from
Ren et al. (2014).
https://academic.oup.com/gji/article/197/2/948/616510
https://library.seg.org/doi/pdf/10.1190/1.2714334
Parameters
----------
p : TimeFunction
Pressure field.
"""
time_order = p.time_order
eq_stencil = stencils[('ren', time_order)]
eqn = eq_stencil(model, geometry, p, forward=forward, **kwargs)
return eqn
def deng_mcmechan(model, geometry, p, forward=True, **kwargs):
"""
Implementation of the 1st order viscoacoustic wave-equation and 2nd order
viscoacoustic wave-equation from Deng and McMechan (2007).
https://library.seg.org/doi/pdf/10.1190/1.2714334
Parameters
----------
p : TimeFunction
Pressure field.
"""
time_order = p.time_order
eq_stencil = stencils[('deng_mcmechan', time_order)]
eqn = eq_stencil(model, geometry, p, forward=forward, **kwargs)
return eqn
def ForwardOperator(model, geometry, space_order=4, kernel='sls', time_order=2,
save=False, **kwargs):
"""
Construct method for the forward modelling operator in a viscoacoustic medium.
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
Space discretization order.
kernel : string, optional
selects a viscoacoustic equation from the options below:
sls (Standard Linear Solid) :
1st order - Blanch and Symes (1995) / Dutta and Schuster (2014)
viscoacoustic equation
2nd order - Bai et al. (2014) viscoacoustic equation
ren - Ren et al. (2014) viscoacoustic equation
deng_mcmechan - Deng and McMechan (2007) viscoacoustic equation
Defaults to sls 2nd order.
save : int or Buffer
Saving flag, True saves all time steps, False saves three buffered
indices (last three time steps). Defaults to False.
"""
# Create symbols for forward wavefield, particle velocity, source and receivers
save_t = geometry.nt if save else None
if time_order == 1:
v = VectorTimeFunction(name="v", grid=model.grid, time_order=time_order,
space_order=space_order, save=save_t)
kwargs.update({'v': v})
p = TimeFunction(name="p", grid=model.grid, time_order=time_order,
space_order=space_order, save=save_t, staggered=NODE)
# Equations kernels
eq_kernel = kernels[kernel]
eqn = eq_kernel(model, geometry, p, save=save, **kwargs)
srcrec = src_rec(p, model, geometry)
# Substitute spacing terms to reduce flops
return Operator(eqn + srcrec, subs=model.spacing_map,
name='Forward', **kwargs)
def AdjointOperator(model, geometry, space_order=4, kernel='sls', time_order=2, **kwargs):
"""
Construct an adjoint modelling operator in a viscoacoustic medium.
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
Space discretization order.
kernel : selects a visco-acoustic equation from the options below:
sls (Standard Linear Solid) :
1st order - Blanch and Symes (1995) / Dutta and Schuster (2014)
viscoacoustic equation
2nd order - Bai et al. (2014) viscoacoustic equation
ren - Ren et al. (2014) viscoacoustic equation
deng_mcmechan - Deng and McMechan (2007) viscoacoustic equation
Defaults to sls 2nd order.
"""
if time_order == 1:
va = VectorTimeFunction(name="va", grid=model.grid, time_order=time_order,
space_order=space_order)
kwargs.update({'v': va})
pa = TimeFunction(name="pa", grid=model.grid, time_order=time_order,
space_order=space_order, staggered=NODE)
# Equations kernels
eq_kernel = kernels[kernel]
eqn = eq_kernel(model, geometry, pa, forward=False, **kwargs)
srcrec = src_rec(pa, model, geometry, forward=False)
# Substitute spacing terms to reduce flops
return Operator(eqn + srcrec, subs=model.spacing_map, name='Adjoint', **kwargs)
def GradientOperator(model, geometry, space_order=4, kernel='sls', time_order=2,
save=True, **kwargs):
"""
Construct a gradient operator in an acoustic media.
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
Space discretization order.
save : int or Buffer, optional
Option to store the entire (unrolled) wavefield.
kernel : selects a visco-acoustic equation from the options below:
sls (Standard Linear Solid) :
1st order - Blanch and Symes (1995) / Dutta and Schuster (2014)
viscoacoustic equation
2nd order - Bai et al. (2014) viscoacoustic equation
ren - Ren et al. (2014) viscoacoustic equation
deng_mcmechan - Deng and McMechan (2007) viscoacoustic equation
Defaults to sls 2nd order.
"""
# Gradient symbol and wavefield symbols
save_t = geometry.nt if save else None
grad = Function(name='grad', grid=model.grid)
p = TimeFunction(name='p', grid=model.grid, time_order=2, space_order=space_order,
save=save_t, staggered=NODE)
pa = TimeFunction(name='pa', grid=model.grid, time_order=time_order,
space_order=space_order, staggered=NODE)
# Equations kernels
eq_kernel = kernels[kernel]
eqn = eq_kernel(model, geometry, pa, forward=False, save=False, **kwargs)
gradient_update = Inc(grad, - p.dt2 * pa)
# Add expression for receiver injection
_, recterm = src_rec(pa, model, geometry, forward=False)
# Substitute spacing terms to reduce flops
return Operator(eqn + recterm + [gradient_update], subs=model.spacing_map,
name='Gradient', **kwargs)
def BornOperator(model, geometry, space_order=4, kernel='sls', time_order=2, **kwargs):
"""
Construct an Linearized Born operator in an acoustic media.
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
Space discretization order.
kernel : str, optional
Type of discretization, centered or shifted.
"""
# Create wavefields and a dm field
p = TimeFunction(name='p', grid=model.grid, time_order=time_order,
space_order=space_order, staggered=NODE)
P = TimeFunction(name='P', grid=model.grid, time_order=time_order,
space_order=space_order, staggered=NODE)
rp = TimeFunction(name="rp", grid=model.grid, time_order=time_order,
space_order=space_order, staggered=NODE)
rP = TimeFunction(name="rP", grid=model.grid, time_order=time_order,
space_order=space_order, staggered=NODE)
dm = Function(name='dm', grid=model.grid, space_order=0)
# Equations kernels
eq_kernel = kernels[kernel]
eqn1 = eq_kernel(model, geometry, p, r=rp, **kwargs)
s = model.grid.stepping_dim.spacing
q = -dm * (p.forward - 2 * p + p.backward) / (s**2)
eqn2 = eq_kernel(model, geometry, P, r=rP, q=q, **kwargs)
# Add source term expression for p
src_term, _ = src_rec(p, model, geometry)
# Create receiver interpolation expression from P
_, rec_term = src_rec(P, model, geometry)
# Substitute spacing terms to reduce flops
return Operator(eqn1 + src_term + rec_term + eqn2, subs=model.spacing_map,
name='Born', **kwargs)
kernels = {'sls': sls, 'ren': ren, 'deng_mcmechan': deng_mcmechan}
stencils = {('sls', 1): sls_1st_order, ('sls', 2): sls_2nd_order,
('deng_mcmechan', 1): deng_1st_order,
('deng_mcmechan', 2): deng_2nd_order,
('ren', 1): ren_1st_order, ('ren', 2): ren_2nd_order}