/
tti_example.py
83 lines (66 loc) · 3.53 KB
/
tti_example.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
import numpy as np
from argparse import ArgumentParser
from devito.logger import warning
from examples.seismic import demo_model, TimeAxis, Receiver, RickerSource
from examples.seismic.tti import AnisotropicWaveSolver
def tti_setup(shape=(50, 50, 50), spacing=(20.0, 20.0, 20.0), tn=250.0,
space_order=4, nbpml=10, **kwargs):
nrec = 101
# Two layer model for true velocity
model = demo_model('layers-tti', shape=shape, spacing=spacing, nbpml=nbpml)
# Derive timestepping from model spacing
dt = model.critical_dt
t0 = 0.0
time_range = TimeAxis(start=t0, stop=tn, step=dt)
# Define source geometry (center of domain, just below surface)
src = RickerSource(name='src', grid=model.grid, f0=0.015, time_range=time_range)
src.coordinates.data[0, :] = np.array(model.domain_size) * .5
src.coordinates.data[0, -1] = model.origin[-1] + 2 * spacing[-1]
# Define receiver geometry (spread across x, lust below surface)
rec = Receiver(name='nrec', grid=model.grid, time_range=time_range, npoint=nrec)
rec.coordinates.data[:, 0] = np.linspace(0., model.domain_size[0], num=nrec)
rec.coordinates.data[:, 1:] = src.coordinates.data[0, 1:]
return AnisotropicWaveSolver(model, source=src, receiver=rec,
space_order=space_order, **kwargs)
def run(shape=(50, 50, 50), spacing=(20.0, 20.0, 20.0), tn=250.0,
autotune=False, time_order=2, space_order=4, nbpml=10,
kernel='centered', **kwargs):
solver = tti_setup(shape, spacing, tn, space_order, nbpml, **kwargs)
if space_order % 4 != 0:
warning('WARNING: TTI requires a space_order that is a multiple of 4!')
rec, u, v, summary = solver.forward(autotune=autotune, kernel=kernel)
return summary.gflopss, summary.oi, summary.timings, [rec, u, v]
if __name__ == "__main__":
description = ("Example script to execute a TTI forward operator.")
parser = ArgumentParser(description=description)
parser.add_argument('--2d', dest='dim2', default=False, action='store_true',
help="Preset to determine the physical problem setup")
parser.add_argument('-a', '--autotune', default=False, action='store_true',
help="Enable autotuning for block sizes")
parser.add_argument("-so", "--space_order", default=4,
type=int, help="Space order of the simulation")
parser.add_argument("--nbpml", default=40,
type=int, help="Number of PML layers around the domain")
parser.add_argument("-k", dest="kernel", default='centered',
choices=['centered', 'shifted'],
help="Choice of finite-difference kernel")
parser.add_argument("-dse", "-dse", default="advanced",
choices=["noop", "basic", "advanced",
"speculative", "aggressive"],
help="Devito symbolic engine (DSE) mode")
parser.add_argument("-dle", default="advanced",
choices=["noop", "advanced", "speculative"],
help="Devito loop engine (DSE) mode")
args = parser.parse_args()
# 3D preset parameters
if args.dim2:
shape = (150, 150)
spacing = (10.0, 10.0)
tn = 750.0
else:
shape = (50, 50, 50)
spacing = (10.0, 10.0, 10.0)
tn = 250.0
run(shape=shape, spacing=spacing, nbpml=args.nbpml, tn=tn,
space_order=args.space_order, autotune=args.autotune, dse=args.dse,
dle=args.dle, kernel=args.kernel)