Skip to content
This repository was archived by the owner on Feb 21, 2023. It is now read-only.

Commit 9f07846

Browse files
committedOct 27, 2021
tests: added 1st-order gather test/plotter
1 parent 1be5c6b commit 9f07846

File tree

1 file changed

+352
-0
lines changed

1 file changed

+352
-0
lines changed
 

‎tests/test_gathers_1st_order.py

+352
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,352 @@
1+
import numpy as np
2+
import pandas as pd
3+
import sympy as sp
4+
import matplotlib.pyplot as plt
5+
import pytest
6+
import sys
7+
8+
from examples.seismic import Model, TimeAxis, RickerSource, Receiver
9+
from devito import TimeFunction, Eq, Operator, NODE, ConditionalDimension, \
10+
Function, Le
11+
from devitoboundary import BoundaryConditions, ImmersedBoundary, AxialDistanceFunction
12+
13+
14+
def reference_shot(model, time_range, f0, s_o):
15+
"""
16+
Produce a reference shot gather with a level, conventional free-surface
17+
implementation.
18+
"""
19+
x, y, z = model.grid.dimensions
20+
t = model.grid.stepping_dim
21+
dt = model.critical_dt/2 # Time step from model grid spacing
22+
23+
vel = 1.5
24+
rho = 1.5
25+
26+
src = RickerSource(name='src', grid=model.grid, f0=f0,
27+
npoint=1, time_range=time_range)
28+
29+
# First, position source centrally in all dimensions, then set depth
30+
src.coordinates.data[0, :] = np.array(model.domain_size) * .5
31+
# Remember that 0, 0, 0 is top left corner
32+
# Depth is 100m from free-surface boundary
33+
src.coordinates.data[0, -1] = 600.
34+
35+
# Create symbol for 101 receivers
36+
rec = Receiver(name='rec', grid=model.grid, npoint=101,
37+
time_range=time_range)
38+
39+
# Prescribe even spacing for receivers along the x-axis
40+
rec.coordinates.data[:, 0] = np.linspace(0, model.domain_size[0], num=101)
41+
rec.coordinates.data[:, 1] = 500. # Centered on y axis
42+
rec.coordinates.data[:, 2] = 650. # Depth is 150m from free surface
43+
44+
# Define the wavefield with the size of the model and the time dimension
45+
p = TimeFunction(name="p", grid=model.grid, time_order=1, space_order=s_o,
46+
staggered=NODE)
47+
vx = TimeFunction(name="v_x", grid=model.grid, time_order=1, space_order=s_o,
48+
staggered=x)
49+
vy = TimeFunction(name="v_y", grid=model.grid, time_order=1, space_order=s_o,
50+
staggered=y)
51+
vz = TimeFunction(name="v_z", grid=model.grid, time_order=1, space_order=s_o,
52+
staggered=z)
53+
54+
# We can now write the system of PDEs
55+
eq_vx = Eq(vx.forward, vx - dt*p.dx/rho)
56+
eq_vy = Eq(vy.forward, vy - dt*p.dy/rho)
57+
eq_vz = Eq(vz.forward, vz - dt*p.dz/rho)
58+
eq_p = Eq(p.forward,
59+
p - dt*vel**2*rho*(vx.forward.dx + vy.forward.dy + vz.forward.dz)/(1+model.damp))
60+
61+
# Finally we define the source injection and receiver read function
62+
src_term = src.inject(field=p.forward,
63+
expr=src)
64+
65+
# Create interpolation expression for receivers
66+
rec_term = rec.interpolate(expr=p.forward)
67+
68+
# Pressure free surface conditions
69+
free_surface_p = [Eq(p[t+1, x, y, 60], 0)]
70+
for i in range(s_o//2):
71+
free_surface_p.append(Eq(p[t+1, x, y, 59-i], -p[t+1, x, y, 61+i]))
72+
73+
free_surface_v = []
74+
for i in range(s_o//2):
75+
free_surface_v.append(Eq(vz[t+1, x, y, 59-i], vz[t+1, x, y, 60+i]))
76+
77+
op = Operator([eq_vx, eq_vy, eq_vz] + free_surface_v
78+
+ [eq_p] + free_surface_p
79+
+ src_term + rec_term)
80+
81+
op(time=time_range.num-1, dt=model.critical_dt/2)
82+
83+
"""
84+
plt.imshow(p.data[-1, :, 60].T, origin='upper')
85+
plt.colorbar()
86+
plt.show()
87+
"""
88+
89+
return rec.data
90+
91+
92+
def tilted_shot(model, time_range, f0, s_o, tilt, qc=False, toggle_normals=False):
93+
"""
94+
Produce a shot for the same setup, but tilted with immersed free surface
95+
"""
96+
# FIXME: Something is wonky here -> am I using a consistent surface?
97+
x, y, z = model.grid.dimensions
98+
h_x, h_y, h_z = model.grid.spacing
99+
dt = model.critical_dt/2 # Time step from model grid spacing
100+
101+
vel = 1.5
102+
rho = 1.
103+
104+
src = RickerSource(name='src', grid=model.grid, f0=f0,
105+
npoint=1, time_range=time_range)
106+
107+
# First, position source, then set depth
108+
src.coordinates.data[0, 0] = 500. - 100.*np.sin(np.radians(tilt))
109+
src.coordinates.data[0, 1] = 500.
110+
# Remember that 0, 0, 0 is top left corner
111+
# Depth is 100m from free-surface boundary
112+
src.coordinates.data[0, 2] = 500. + 100.*np.cos(np.radians(tilt))
113+
114+
# Create symbol for 101 receivers
115+
rec = Receiver(name='rec', grid=model.grid, npoint=101,
116+
time_range=time_range)
117+
118+
# Prescribe even spacing for receivers along the x-axis
119+
rec_center_x = 500. - 150.*np.sin(np.radians(tilt))
120+
rec_center_z = 500. + 150.*np.cos(np.radians(tilt))
121+
122+
rec_top_x = rec_center_x - 500.*np.cos(np.radians(tilt))
123+
rec_bottom_x = rec_center_x + 500.*np.cos(np.radians(tilt))
124+
125+
rec_top_z = rec_center_z - 500.*np.sin(np.radians(tilt))
126+
rec_bottom_z = rec_center_z + 500.*np.sin(np.radians(tilt))
127+
128+
rec.coordinates.data[:, 0] = np.linspace(rec_top_x, rec_bottom_x, num=101)
129+
rec.coordinates.data[:, 1] = 500. # Centered on y axis
130+
rec.coordinates.data[:, 2] = np.linspace(rec_top_z, rec_bottom_z, num=101)
131+
132+
# Define the wavefield with the size of the model and the time dimension
133+
p = TimeFunction(name="p", grid=model.grid, time_order=2, space_order=s_o,
134+
staggered=NODE, coefficients='symbolic')
135+
vx = TimeFunction(name="v_x", grid=model.grid, time_order=2, space_order=s_o+2,
136+
staggered=x, coefficients='symbolic')
137+
vy = TimeFunction(name="v_y", grid=model.grid, time_order=2, space_order=s_o+2,
138+
staggered=y, coefficients='symbolic')
139+
vz = TimeFunction(name="v_z", grid=model.grid, time_order=2, space_order=s_o+2,
140+
staggered=z, coefficients='symbolic')
141+
142+
# Dummy functions
143+
vx_d = TimeFunction(name='v_x_d', grid=model.grid,
144+
space_order=s_o, staggered=x, coefficients='symbolic')
145+
vy_d = TimeFunction(name='v_y_d', grid=model.grid,
146+
space_order=s_o, staggered=y, coefficients='symbolic')
147+
vz_d = TimeFunction(name='v_z_d', grid=model.grid,
148+
space_order=s_o, staggered=z, coefficients='symbolic')
149+
150+
infile = 'tests/trial_surfaces/angled_surface_'+str(tilt)+'.ply'
151+
152+
# Zero even pressure derivatives on the boundary
153+
spec_p = {2*i: 0 for i in range(p.space_order)}
154+
bcs_p = BoundaryConditions(spec_p, p.space_order)
155+
# Zero odd derivatives on the boundary
156+
spec_v = {2*i+1: 0 for i in range(vx_d.space_order)}
157+
bcs_v = BoundaryConditions(spec_v, vx_d.space_order)
158+
159+
functions = pd.DataFrame({'function': [p, vx_d, vy_d, vz_d],
160+
'bcs': [bcs_p, bcs_v, bcs_v, bcs_v],
161+
'subs_function': [None, vx, vy, vz]},
162+
columns=['function', 'bcs', 'subs_function'])
163+
164+
# Create the immersed boundary surface
165+
surface = ImmersedBoundary('topography', infile, functions,
166+
interior_point=tuple(src.coordinates.data[0]),
167+
qc=qc, toggle_normals=toggle_normals)
168+
169+
# Configure derivative needed
170+
derivs = pd.DataFrame({'function': [p, vx_d, vy_d, vz_d],
171+
'derivative': [1, 1, 1, 1],
172+
'eval_offset': [(0.5, 0.5, 0.5), (-0.5, 0.5, 0.5),
173+
(0.5, -0.5, 0.5), (0.5, 0.5, -0.5)]},
174+
columns=['function', 'derivative', 'eval_offset'])
175+
176+
more_derivs = pd.DataFrame({'function': [p],
177+
'derivative': [2],
178+
'eval_offset': [(0., 0., 0.)]},
179+
columns=['function', 'derivative',
180+
'eval_offset'])
181+
182+
coeffs = surface.subs(derivs)
183+
more_coeffs = surface.subs(more_derivs)
184+
185+
# Set up AxialDistanceFunction and conditional dimensions so that the
186+
# 2nd-order pressure update is applied when |eta| < 0.5
187+
ax = AxialDistanceFunction(p, infile,
188+
toggle_normals=toggle_normals)
189+
190+
# Needed to sidestep a compilation-preventing bug in Devito
191+
ax_x = Function(name='ax_x', shape=model.grid.shape, dimensions=(x, y, z))
192+
ax_y = Function(name='ax_y', shape=model.grid.shape, dimensions=(x, y, z))
193+
ax_z = Function(name='ax_z', shape=model.grid.shape, dimensions=(x, y, z))
194+
195+
ax_x.data[:] = np.array(ax.axial[0].data[:])
196+
ax_y.data[:] = np.array(ax.axial[1].data[:])
197+
ax_z.data[:] = np.array(ax.axial[2].data[:])
198+
199+
second_update = sp.Or(sp.Or(Le(sp.Abs(ax_x), 0.5*h_x),
200+
Le(sp.Abs(ax_y), 0.5*h_y)),
201+
Le(sp.Abs(ax_z), 0.5*h_z))
202+
203+
# Conditional masks for update
204+
use_2nd = ConditionalDimension(name='use_2nd', parent=z,
205+
condition=second_update)
206+
use_1st = ConditionalDimension(name='use_1st', parent=z,
207+
condition=sp.Not(second_update))
208+
209+
# We can now write the system of PDEs
210+
eq_vx = Eq(vx.forward, vx - dt*p.dx/rho, coefficients=coeffs)
211+
eq_vy = Eq(vy.forward, vy - dt*p.dy/rho, coefficients=coeffs)
212+
eq_vz = Eq(vz.forward, vz - dt*p.dz/rho, coefficients=coeffs)
213+
eq_p = Eq(p.forward,
214+
p - dt*vel**2*rho*(vx.forward.dx + vy.forward.dy + vz.forward.dz)/(1+model.damp), coefficients=coeffs, implicit_dims=use_1st)
215+
eq_p2 = Eq(p.forward,
216+
2*p - p.backward + dt**2*vel**2*p.laplace/(1+model.damp),
217+
coefficients=more_coeffs, implicit_dims=use_2nd)
218+
219+
# Finally we define the source injection and receiver read function
220+
src_term = src.inject(field=p.forward,
221+
expr=src)
222+
223+
# Create interpolation expression for receivers
224+
rec_term = rec.interpolate(expr=p.forward)
225+
226+
op = Operator([eq_vx, eq_vy, eq_vz, eq_p, eq_p2] + src_term + rec_term)
227+
op(time=time_range.num-1, dt=model.critical_dt/2)
228+
229+
"""
230+
plt.imshow(p.data[-1, :, 60].T, origin='upper')
231+
plt.colorbar()
232+
plt.show()
233+
"""
234+
235+
return rec.data
236+
237+
238+
class TestGathers:
239+
"""
240+
A class for testing the accuracy of gathers resulting from a reflection off
241+
the immersed boundary.
242+
"""
243+
244+
@pytest.mark.parametrize('s_o', [4, 6])
245+
@pytest.mark.parametrize('spec', [(5, False), (10, True), (15, True),
246+
(20, True), (25, True), (30, True),
247+
(35, True), (40, True), (45, False)])
248+
def test_tilted_boundary(self, s_o, spec):
249+
"""
250+
Check that gathers for a tilted boundary match those generated with a
251+
conventional horizontal free surface and the same geometry.
252+
"""
253+
tilt, toggle_normals = spec
254+
max_thres = 0.0026
255+
avg_thres = 2.5e-4
256+
257+
# Define a physical size
258+
shape = (101, 101, 101) # Number of grid point (nx, ny, nz)
259+
spacing = (10., 10., 10.) # Grid spacing in m. The domain size is 1x1x1km
260+
origin = (0., 0., 0.)
261+
262+
v = 1.5
263+
264+
model = Model(vp=v, origin=origin, shape=shape, spacing=spacing,
265+
space_order=4, nbl=10, bcs="damp")
266+
267+
t0 = 0. # Simulation starts a t=0
268+
tn = 500. # Simulation last 0.5 seconds (500 ms)
269+
dt = model.critical_dt/2 # Time step from model grid spacing
270+
271+
time_range = TimeAxis(start=t0, stop=tn, step=dt)
272+
273+
f0 = 0.010 # Source peak frequency is 10Hz (0.010 kHz)
274+
275+
ref = reference_shot(model, time_range, f0, s_o)
276+
tilted = tilted_shot(model, time_range, f0, s_o, tilt,
277+
toggle_normals=toggle_normals, qc=False)
278+
279+
"""
280+
plt.imshow(ref, aspect='auto', vmin=-0.03, vmax=0.03, cmap='seismic')
281+
plt.colorbar()
282+
plt.show()
283+
284+
plt.imshow(tilted, aspect='auto', vmin=-0.03, vmax=0.03, cmap='seismic')
285+
plt.colorbar()
286+
plt.show()
287+
288+
plt.imshow(tilted-ref, aspect='auto', vmin=-0.03, vmax=0.03, cmap='seismic')
289+
plt.colorbar()
290+
plt.show()
291+
"""
292+
293+
assert np.amax(np.absolute(ref - tilted)) < max_thres
294+
assert np.mean(np.absolute(ref-tilted)) < avg_thres
295+
296+
297+
def main(s_o, tilt, filepath, toggle_normals=True, qc=True):
298+
# Define a physical size
299+
shape = (101, 101, 101) # Number of grid point (nx, ny, nz)
300+
spacing = (10., 10., 10.) # Grid spacing in m. The domain size is 1x1x1km
301+
origin = (0., 0., 0.)
302+
303+
v = 1.5
304+
305+
model = Model(vp=v, origin=origin, shape=shape, spacing=spacing,
306+
space_order=4, nbl=10, bcs="damp")
307+
308+
t0 = 0. # Simulation starts a t=0
309+
tn = 500. # Simulation last 0.5 seconds (500 ms)
310+
dt = model.critical_dt/2 # Time step from model grid spacing
311+
312+
time_range = TimeAxis(start=t0, stop=tn, step=dt)
313+
314+
f0 = 0.010 # Source peak frequency is 10Hz (0.010 kHz)
315+
316+
ref = reference_shot(model, time_range, f0, s_o)
317+
tilted = tilted_shot(model, time_range, f0, s_o, tilt,
318+
toggle_normals=toggle_normals, qc=qc)
319+
320+
fig, (ax0, ax1, ax2) = plt.subplots(nrows=1, ncols=3, figsize=(20, 6))
321+
322+
im0 = ax0.imshow(ref, aspect='auto', cmap='seismic', extent=[-500, 500, tn, t0],
323+
vmin=-0.03, vmax=0.03) # noqa: F841
324+
ax0.set_title("Reference")
325+
ax0.set_xlabel("Offset (m)")
326+
ax0.set_ylabel("Time (ms)")
327+
328+
im1 = ax1.imshow(tilted, aspect='auto', cmap='seismic', extent=[-500, 500, tn, t0],
329+
vmin=-0.03, vmax=0.03) # noqa: F841
330+
ax1.set_title("{} degree tilt".format(str(tilt)))
331+
ax1.set_xlabel("Offset (m)")
332+
333+
im2 = ax2.imshow(ref - tilted, aspect='auto', cmap='seismic', extent=[-500, 500, tn, t0],
334+
vmin=-0.03, vmax=0.03) # noqa: F841
335+
ax2.set_title("Difference")
336+
ax2.set_xlabel("Offset (m)")
337+
fig.colorbar(im2)
338+
fig.tight_layout()
339+
if filepath == 'show':
340+
plt.show()
341+
else:
342+
plt.savefig(filepath + "/order_{}_tilt_{}".format(str(s_o), str(tilt)))
343+
plt.close()
344+
345+
346+
if __name__ == "__main__":
347+
s_o = int(sys.argv[1])
348+
tilt = int(sys.argv[2])
349+
filepath = sys.argv[3]
350+
toggle_normals = bool(int(sys.argv[4]))
351+
qc = bool(int(sys.argv[5]))
352+
main(s_o, tilt, filepath, toggle_normals=toggle_normals, qc=qc)

0 commit comments

Comments
 (0)
This repository has been archived.