-
Notifications
You must be signed in to change notification settings - Fork 0
/
observationoperator.py
369 lines (315 loc) · 12.5 KB
/
observationoperator.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
"""
Define observation operators for inverse problem
"""
import sys
import abc
import numpy as np
from numpy import sqrt
#from numpy.linalg import norm
#from numpy.random import randn
try:
import matplotlib.pyplot as plt
except:
pass
from dolfin import Function, TrialFunction, TestFunction, \
Constant, Point, as_backend_type, assemble, inner, dx
from miscfenics import isFunction, isVector, isarray, arearrays, isequal, setfct
from sourceterms import PointSources
class TimeFilter():
""" Create time filter to fade out data misfit (hence src term in adj eqn) """
def __init__(self, times=None):
""" Input times:
times[0] = t0 = initial time
times[1] = t1 = beginning of flat section at 1.0
times[2] = t2 = end of flat section at 1.0
times[3] = T = final time """
if times == None:
self.t0 = -np.inf
self.t1 = -np.inf
self.t2 = np.inf
self.T = np.inf
else:
self.t0 = times[0]
self.t1 = times[1]
self.t2 = times[2]
self.T = times[3]
self.t1b = 0.5*(self.t0 + self.t1)
self.t2b = 0.5*(self.t2 + self.T)
#self.t2b = 2*self.t2-self.T
def __call__(self, tt):
""" Overload () operator """
#assert tt >= self.t0 and tt <= self.T, \
#'Input tt={} out of bounds [t0, T]; tt-t0={}, T-tt={}'.format(tt, tt-self.t0, self.T-tt)
if tt <= self.t0 + 1e-16: return 0.0
if tt >= self.T - 1e-16: return 0.0
if tt <= self.t1b:
return 0.5*np.exp(1./(self.t0-tt))/np.exp(1./(self.t0-self.t1b))
if tt <= self.t1:
return 1.0 - self.__call__(tt-2.0*(tt-self.t1b))
if tt >= self.t2b:
return 0.5*np.exp(1./(tt-self.T))/np.exp(1./(self.t2b-self.T))
if tt >= self.t2:
return 1.0 - self.__call__(self.T-(tt-self.t2))
return 1.0
def evaluate(self, times):
""" vectorized verstion of __call__ """
vectcall = np.vectorize(self.__call__)
return vectcall(times)
def plot(self, ndt=1000):
""" Plot the shape of the filter along with its fft """
tt = np.linspace(self.t0, self.T, ndt)
xx = self.evaluate(tt)
ff = np.fft.fft(xx) # fft(time-domain)
ffn = np.sqrt(ff.real**2 + ff.imag**2)
ffxi = np.fft.fftfreq(len(xx), d=tt[1]-tt[0])
fig = plt.figure()
ax1 = fig.add_subplot(121)
ax1.plot(tt, xx)
ax2 = fig.add_subplot(122)
ax2.plot(np.fft.fftshift(ffxi), np.fft.fftshift(ffn))
return fig
###########################################################
###########################################################
class ObservationOperator():
"""Define observation operator and all actions using this observation
operator, i.e., cost function, rhs in adj eqn and term in incr. adj eqn"""
__metaclass__ = abc.ABCMeta
#Instantiation
def __init__(self, parameters, mycomm):
self.parameters = parameters
if self.parameters.has_key('noise'):
self.noise = True
self.noisepercent = self.parameters['noise']
else:
self.noise = False
self.mycomm = mycomm
self._assemble()
@abc.abstractmethod
def _assemble(self): print "Needs to be implemented"
@abc.abstractmethod
def obs(self, uin): print "Needs to be implemented"
# Note: This must return the global value (in parallel)
@abc.abstractmethod
def costfct(self, uin, udin): print "Needs to be implemented"
@abc.abstractmethod
def assemble_rhsadj(self, uin, udin): print "Needs to be implemented"
@abc.abstractmethod
def incradj(self, uin): print "Needs to be implemented"
# def apply_noise(self, uin):
# """Apply Gaussian noise to np.array of data.
# noisepercent = 0.02 => 2% noise level, i.e.,
# || u - ud || / || ud || = || noise || / || ud || = 0.02"""
# isarray(uin)
# noisevect = randn(len(uin))
# # Get norm of entire random vector and
# # Get norm of entire vector ud (not just local part):
# normrand = sqrt(MPI.sum(self.mycomm, norm(noisevect)**2))
# normud = sqrt(MPI.sum(self.mycomm, norm(uin)**2))
# noisevect /= normrand
# noisevect *= self.noisepercent * normud
# objnoise_glob = (self.noisepercent * normud)**2
# UDnoise = uin + noisevect
# return UDnoise, objnoise_glob
###########################################################
# Derived Classes
###########################################################
class ObsEntireDomain(ObservationOperator):
"""Observation operator over the entire domain
with L2 misfit norm
parameters must be dictionary containing:
V = function space for state variable"""
def _assemble(self):
self.V = self.parameters['V']
self.diff = Function(self.V)
self.diffv = self.diff.vector()
self.trial = TrialFunction(self.V)
self.test = TestFunction(self.V)
self.W = assemble(inner(self.trial, self.test)*dx)
def obs(self, uin):
isFunction(uin)
if not(self.noise): return uin.vector().array(), 0.0
else: return self.apply_noise(uin.vector().array())
# Note: this returns global value (in parallel)
def costfct(self, uin, udin):
arearrays(uin, udin)
self.diffv[:] = uin - udin
return 0.5 * (self.W*self.diffv).inner(self.diffv)
def assemble_rhsadj(self, uin, udin, outp, bc):
arearrays(uin, udin)
isFunction(outp)
self.diffv[:] = uin - udin
outp.vector()[:] = - (self.W * self.diffv).array()
if bc is not None:
bc.apply(outp.vector())
def incradj(self, uin):
isFunction(uin)
return self.hessian(uin.vector())
# Make non-array version of cost:
def costfct_F(self, uin, udin):
isFunction(uin)
isFunction(udin)
setfct(self.diff, uin.vector()-udin.vector())
return 0.5 * (self.W*self.diffv).inner(self.diffv)
def grad(self, uin, udin):
isFunction(uin)
isFunction(udin)
setfct(self.diff, uin.vector() - udin.vector())
return self.W * self.diffv
def hessian(self, uin):
isVector(uin)
return self.W * uin
##################
class ObsPointwise(ObservationOperator):
"""Observation operator at finite nb of points
parameters must be a dictionary containing:
V = function space for state variable
Points = list of coordinates"""
def _assemble(self):
self.V = self.parameters['V']
self.Points = self.parameters['Points']
self.nbPts = len(self.Points)
self.BtBu = Function(self.V)
PtSrc = PointSources(self.V, self.Points)
self.B = PtSrc.PtSrc
def Bdotlocal(self, uin):
"""Compute B.uin as a np.array, using only local info
uin must be a Function(self.V)"""
isFunction(uin)
Bu = np.zeros(self.nbPts)
for ii, bb in enumerate(self.B):
Bu[ii] = np.dot(bb.array(), uin.vector().array()) # Note: local inner-product
return Bu
def Bdot(self, uin):
"""Compute B.uin as a np.array, using global info
uin must be a Function(self.V)"""
isFunction(uin)
Bu = np.zeros(self.nbPts)
for ii, bb in enumerate(self.B):
Bu[ii] = bb.inner(uin.vector()) # Note: global inner-product
return Bu
def BTdot(self, uin):
"""Compute B^T.uin as a np.array
uin must be a np.array"""
isarray(uin)
u = Function(self.V)
out = u.vector()
for ii, bb in enumerate(self.B):
out.axpy(uin[ii], bb)
return out.array()
#@profile
def BTdotvec(self, uin, outvect):
""" Compute B^T.uin """
isarray(uin)
outvect.zero()
for ii, bb in enumerate(self.B):
outvect.axpy(uin[ii], bb)
def obs(self, uin):
"""Compute B.uin + eps, where eps is noise
uin must be a Function(V)"""
if not(self.noise): return self.Bdot(uin), 0.0
else:
sys.exit(1)
# Bref = self.Bdot(uin)
# uin_noise, tmp = self.apply_noise(uin.vector().array())
# unoise = Function(self.V)
# unoise.vector()[:] = uin_noise
# Bnoise = self.Bdot(unoise)
# diff = Bref - Bnoise
# noiselevel = np.dot(diff, diff)
# try:
# noiselevel_glob = MPI.sum(self.mycomm, noiselevel)
# except:
# noiselevel_glob = noiselevel
# return Bnoise, noiselevel_glob
def costfct(self, uin, udin):
"""Compute cost functional from observed fwd and data, i.e.,
return .5*||uin - udin||^2.
uin & udin are np.arrays"""
arearrays(uin, udin)
diff = uin - udin
return 0.5*np.dot(diff, diff)
def assemble_rhsadj(self, uin, udin, outp, bc):
"""Compute rhs term for adjoint equation and store it in outp, i.e.,
outp = - B^T( uin - udin), where uin = obs(fwd solution)
uin & udin = np.arrays
outp = Function(self.V)
bc = fenics' boundary conditons"""
arearrays(uin, udin)
isFunction(outp)
diff = uin - udin
outp.vector()[:] = -1.0 * self.BTdot(diff)
if bc is not None:
bc.apply(outp.vector())
def incradj(self, uin):
"""Compute the observation part of the incremental adjoint equation, i.e,
return B^T.B.uin
uin = Function(self.V)"""
isFunction(uin)
self.BtBu.vector()[:] = self.BTdot( self.Bdotlocal(uin) )
return self.BtBu.vector()
##################
class TimeObsPtwise():
"""
Create time-dependent pointwise observation operator
Arguments to assemble:
paramObsPtwise = parameters used to create ptwise observation operator
timefilter = [t0, t1, t2, T], times to initialize time-filtering function
"""
def __init__(self, paramObsPtwise, timefilter=None):
mycomm = paramObsPtwise['V'].mesh().mpi_comm()
self.PtwiseObs = ObsPointwise(paramObsPtwise, mycomm)
u = Function(self.PtwiseObs.V)
self.outvec = u.vector()
self.st = TimeFilter(timefilter)
def obs(self, uin):
""" return result from pointwise observation w/o time-filtering """
return self.PtwiseObs.Bdot(uin)
def costfct(self, uin, udin, times, factors):
""" compute cost functional
int_0^T s(t) | B u - udin |^2 dt
uin, udin must be in np.array format of shape 'recv x time'
times should be array containg values of t0, t1,..., T,
and time-steps should be all equal (!!)
"""
assert uin.shape == udin.shape, "uin and udin must have same shape"
assert uin.shape[0] == len(times) or uin.shape[1] == len(times), \
"must have as many time steps in uin as in times"
if uin.shape[0] == len(times): diffuinudinsq = (uin.T - udin.T)**2
else: diffuinudinsq = (uin - udin)**2
#
diff = diffuinudinsq*factors*(self.st.evaluate(times))
return 0.5*(diff.sum().sum())
def assemble_rhsadj(self, uin, udin, times, bcadj):
""" Assemble data for rhs of adj eqn
uin, udin = observations and data @ receivers
times = time steps of solution
bcadj = boundary conditions
uin.shape must be 'recv x time' """
assert uin.shape == udin.shape, "uin and udin must have same shape"
assert uin.shape[0] == len(times) or uin.shape[1] == len(times), \
"must have as many time steps in uin as in times"
#
if uin.shape[0] == len(times):
self.diff = uin.T - udin.T
else:
self.diff = uin - udin
self.times = times
self.bcadj = bcadj
#@profile
def ftimeadj(self, tt):
""" Evaluate source term for adj eqn at time tt """
try:
index = int(np.where(isequal(self.times, tt, 1e-14))[0])
except:
print 'Error in ftimeadj at time {}'.format(tt)
print np.min(np.abs(self.times-tt))
sys.exit(0)
dd = self.diff[:, index]
self.PtwiseObs.BTdotvec(dd, self.outvec)
if not self.bcadj == None: self.bcadj.apply(self.outvec.vector())
return self.outvec*(-1.0*self.st(tt))
#@profile
def incradj(self, uhat, tt):
""" Compute B^T B uhat """
self.PtwiseObs.BTdotvec(self.obs(uhat), self.outvec)
return self.outvec*self.st(tt)