-
Notifications
You must be signed in to change notification settings - Fork 221
/
test_gradient.py
212 lines (179 loc) · 9.04 KB
/
test_gradient.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
import numpy as np
import pytest
from numpy import linalg
from conftest import skipif
from devito import Function, info, clear_cache
from examples.seismic.acoustic.acoustic_example import smooth, acoustic_setup as setup
from examples.seismic import Receiver
pytestmark = skipif(['yask', 'ops'])
class TestGradient(object):
def setup_method(self, method):
# Some of these tests are memory intensive as it requires to store the entire
# forward wavefield to compute the gradient (nx.ny.nz.nt). We therefore call
# 'clear_cache()' to release any remaining memory from the previous tests or
# previous instances (different parametrizations) of these tests
clear_cache()
@pytest.mark.parametrize('space_order', [4])
@pytest.mark.parametrize('kernel', ['OT2'])
@pytest.mark.parametrize('shape', [(70, 80)])
def test_gradient_checkpointing(self, shape, kernel, space_order):
"""
This test ensures that the FWI gradient computed with devito
satisfies the Taylor expansion property:
.. math::
\Phi(m0 + h dm) = \Phi(m0) + \O(h) \\
\Phi(m0 + h dm) = \Phi(m0) + h \nabla \Phi(m0) + \O(h^2) \\
\Phi(m0) = .5* || F(m0 + h dm) - D ||_2^2
where
.. math::
\nabla \Phi(m0) = <J^T \delta d, dm> \\
\delta d = F(m0+ h dm) - D \\
with F the Forward modelling operator.
:param dimensions: size of the domain in all dimensions
in number of grid points
:param time_order: order of the time discretization scheme
:param space_order: order of the spacial discretization scheme
:return: assertion that the Taylor properties are satisfied
"""
spacing = tuple(10. for _ in shape)
wave = setup(shape=shape, spacing=spacing, dtype=np.float64,
kernel=kernel, space_order=space_order,
nbpml=40)
m0 = Function(name='m0', grid=wave.model.grid, space_order=space_order)
smooth(m0, wave.model.m)
# Compute receiver data for the true velocity
rec, u, _ = wave.forward()
# Compute receiver data and full wavefield for the smooth velocity
rec0, u0, _ = wave.forward(m=m0, save=True)
# Gradient: <J^T \delta d, dm>
residual = Receiver(name='rec', grid=wave.model.grid, data=rec0.data - rec.data,
time_range=wave.geometry.time_axis,
coordinates=wave.geometry.rec_positions)
gradient, _ = wave.gradient(residual, u0, m=m0, checkpointing=True)
gradient2, _ = wave.gradient(residual, u0, m=m0, checkpointing=False)
assert np.allclose(gradient.data, gradient2.data)
@pytest.mark.parametrize('space_order', [4])
@pytest.mark.parametrize('kernel', ['OT2'])
@pytest.mark.parametrize('shape', [(70, 80)])
@pytest.mark.parametrize('checkpointing', [True, False])
def test_gradientFWI(self, shape, kernel, space_order, checkpointing):
"""
This test ensures that the FWI gradient computed with devito
satisfies the Taylor expansion property:
.. math::
\Phi(m0 + h dm) = \Phi(m0) + \O(h) \\
\Phi(m0 + h dm) = \Phi(m0) + h \nabla \Phi(m0) + \O(h^2) \\
\Phi(m0) = .5* || F(m0 + h dm) - D ||_2^2
where
.. math::
\nabla \Phi(m0) = <J^T \delta d, dm> \\
\delta d = F(m0+ h dm) - D \\
with F the Forward modelling operator.
:param dimensions: size of the domain in all dimensions
in number of grid points
:param time_order: order of the time discretization scheme
:param space_order: order of the spacial discretization scheme
:return: assertion that the Taylor properties are satisfied
"""
spacing = tuple(10. for _ in shape)
wave = setup(shape=shape, spacing=spacing, dtype=np.float64,
kernel=kernel, space_order=space_order,
nbpml=40)
m0 = Function(name='m0', grid=wave.model.grid, space_order=space_order)
smooth(m0, wave.model.m)
dm = np.float32(wave.model.m.data[:] - m0.data[:])
# Compute receiver data for the true velocity
rec, u, _ = wave.forward()
# Compute receiver data and full wavefield for the smooth velocity
rec0, u0, _ = wave.forward(m=m0, save=True)
# Objective function value
F0 = .5*linalg.norm(rec0.data - rec.data)**2
# Gradient: <J^T \delta d, dm>
residual = Receiver(name='rec', grid=wave.model.grid, data=rec0.data - rec.data,
time_range=wave.geometry.time_axis,
coordinates=wave.geometry.rec_positions)
gradient, _ = wave.gradient(residual, u0, m=m0, checkpointing=checkpointing)
G = np.dot(gradient.data.reshape(-1), dm.reshape(-1))
# FWI Gradient test
H = [0.5, 0.25, .125, 0.0625, 0.0312, 0.015625, 0.0078125]
error1 = np.zeros(7)
error2 = np.zeros(7)
for i in range(0, 7):
# Add the perturbation to the model
def initializer(data):
data[:] = m0.data + H[i] * dm
mloc = Function(name='mloc', grid=wave.model.m.grid, space_order=space_order,
initializer=initializer)
# Data for the new model
d = wave.forward(m=mloc)[0]
# First order error Phi(m0+dm) - Phi(m0)
F_i = .5*linalg.norm((d.data - rec.data).reshape(-1))**2
error1[i] = np.absolute(F_i - F0)
# Second order term r Phi(m0+dm) - Phi(m0) - <J(m0)^T \delta d, dm>
error2[i] = np.absolute(F_i - F0 - H[i] * G)
# Test slope of the tests
p1 = np.polyfit(np.log10(H), np.log10(error1), 1)
p2 = np.polyfit(np.log10(H), np.log10(error2), 1)
info('1st order error, Phi(m0+dm)-Phi(m0): %s' % (p1))
info('2nd order error, Phi(m0+dm)-Phi(m0) - <J(m0)^T \delta d, dm>: %s' % (p2))
assert np.isclose(p1[0], 1.0, rtol=0.1)
assert np.isclose(p2[0], 2.0, rtol=0.1)
@pytest.mark.parametrize('space_order', [4])
@pytest.mark.parametrize('kernel', ['OT2'])
@pytest.mark.parametrize('shape', [(70, 80)])
def test_gradientJ(self, shape, kernel, space_order):
"""
This test ensures that the Jacobian computed with devito
satisfies the Taylor expansion property:
.. math::
F(m0 + h dm) = F(m0) + \O(h) \\
F(m0 + h dm) = F(m0) + J dm + \O(h^2) \\
with F the Forward modelling operator.
:param dimensions: size of the domain in all dimensions
in number of grid points
:param time_order: order of the time discretization scheme
:param space_order: order of the spacial discretization scheme
:return: assertion that the Taylor properties are satisfied
"""
spacing = tuple(15. for _ in shape)
wave = setup(shape=shape, spacing=spacing, dtype=np.float64,
kernel=kernel, space_order=space_order,
tn=1000., nbpml=10+space_order/2)
m0 = Function(name='m0', grid=wave.model.grid, space_order=space_order)
smooth(m0, wave.model.m)
dm = np.float64(wave.model.m.data - m0.data)
linrec = Receiver(name='rec', grid=wave.model.grid,
time_range=wave.geometry.time_axis,
coordinates=wave.geometry.rec_positions)
# Compute receiver data and full wavefield for the smooth velocity
rec, u0, _ = wave.forward(m=m0, save=False)
# Gradient: J dm
Jdm, _, _, _ = wave.born(dm, rec=linrec, m=m0)
# FWI Gradient test
H = [0.5, 0.25, .125, 0.0625, 0.0312, 0.015625, 0.0078125]
error1 = np.zeros(7)
error2 = np.zeros(7)
for i in range(0, 7):
# Add the perturbation to the model
def initializer(data):
data[:] = m0.data + H[i] * dm
mloc = Function(name='mloc', grid=wave.model.m.grid, space_order=space_order,
initializer=initializer)
# Data for the new model
d = wave.forward(m=mloc)[0]
delta_d = (d.data - rec.data).reshape(-1)
# First order error F(m0 + hdm) - F(m0)
error1[i] = np.linalg.norm(delta_d, 1)
# Second order term F(m0 + hdm) - F(m0) - J dm
error2[i] = np.linalg.norm(delta_d - H[i] * Jdm.data.reshape(-1), 1)
# Test slope of the tests
p1 = np.polyfit(np.log10(H), np.log10(error1), 1)
p2 = np.polyfit(np.log10(H), np.log10(error2), 1)
info('1st order error, Phi(m0+dm)-Phi(m0) with slope: %s compared to 1' % (p1[0]))
info('2nd order error, Phi(m0+dm)-Phi(m0) - <J(m0)^T \delta d, dm>with slope:'
' %s comapred to 2' % (p2[0]))
assert np.isclose(p1[0], 1.0, rtol=0.1)
assert np.isclose(p2[0], 2.0, rtol=0.1)
if __name__ == "__main__":
TestGradient().test_gradientFWI(shape=(70, 80), kernel='OT2', space_order=4,
checkpointing=False)