-
Notifications
You must be signed in to change notification settings - Fork 4
/
__init__.py
421 lines (314 loc) · 14 KB
/
__init__.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
"""Adams-Bashforth ODE solvers."""
__copyright__ = """
Copyright (C) 2007 Andreas Kloeckner
Copyright (C) 2014, 2015 Matt Wala
Copyright (C) 2015 Cory Mikida
"""
__license__ = """
Permission is hereby granted, free of charge, to any person obtaining a copy
of this software and associated documentation files (the "Software"), to deal
in the Software without restriction, including without limitation the rights
to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
copies of the Software, and to permit persons to whom the Software is
furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in
all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
THE SOFTWARE.
"""
import numpy as np
import numpy.linalg as la
from pymbolic import var
from leap import MethodBuilder
__doc__ = """
.. autoclass:: AdamsIntegrationFunctionFamily
.. autoclass:: AdamsMonomialIntegrationFunctionFamily
.. autoclass:: AdamsBashforthMethodBuilder
"""
# {{{ Adams-Bashforth integration (with and without dynamic time steps)
def _linear_comb(coefficients, vectors):
from functools import reduce
from operator import add
return reduce(add,
(coeff * v for coeff, v in zip(coefficients, vectors)))
class AdamsIntegrationFunctionFamily:
"""An abstract interface for function families used for
Adams-type time integration.
.. automethod:: __len__
.. automethod:: evaluate
.. automethod:: antiderivative
"""
def __len__(self):
raise NotImplementedError()
def evaluate(self, func_idx, x):
raise NotImplementedError()
def antiderivative(self, func_idx, x):
raise NotImplementedError()
class AdamsMonomialIntegrationFunctionFamily(AdamsIntegrationFunctionFamily):
"""
Implements :class:`AdamsMonomialIntegrationFunctionFamily`.
"""
def __init__(self, order):
self.order = order
def __len__(self):
return self.order
def evaluate(self, func_idx, x):
return x**func_idx
def antiderivative(self, func_idx, x):
return 1/(func_idx+1) * x**(func_idx+1)
def _emit_func_family_operation(cb, name_gen,
function_family, time_values, hist_vars, rhs_func):
if isinstance(time_values, var):
# {{{ variable time step
hist_len = len(hist_vars)
nfunctions = len(function_family)
array = var("<builtin>array")
linear_solve = var("<builtin>linear_solve")
svd = var("<builtin>svd")
matmul = var("<builtin>matmul")
transpose = var("<builtin>transpose")
# use:
# Vandermonde^T * a_coeffs = integrate(t_start, t_end, monomials)
vdmt = var(name_gen("vdm_transpose"))
cb(vdmt, array(nfunctions*hist_len))
coeff_rhs = var(name_gen("coeff_rhs"))
cb(coeff_rhs, array(nfunctions))
j = var(name_gen("vdm_j"))
for i in range(len(function_family)):
cb(vdmt[i + j*nfunctions], function_family.evaluate(i, time_values[j]),
loops=[(j.name, 0, hist_len)])
for i in range(len(function_family)):
cb(coeff_rhs[i], rhs_func(i))
a_coeffs = var(name_gen("a_coeffs"))
if hist_len == nfunctions:
cb(a_coeffs, linear_solve(vdmt, coeff_rhs, nfunctions, 1))
else:
# Least squares with SVD builtin
u = var(name_gen("u"))
ut = var(name_gen("ut"))
intermed = var(name_gen("intermed"))
ainv = var(name_gen("ainv"))
sigma = var(name_gen("sigma"))
sig_array = var(name_gen("sig_array"))
v = var(name_gen("v"))
vt = var(name_gen("vt"))
cb(ainv, array(nfunctions*hist_len))
cb(intermed, array(nfunctions*hist_len))
cb((u, sigma, vt), svd(vdmt, hist_len))
cb(ut, transpose(u, nfunctions))
cb(v, transpose(vt, hist_len))
# Make singular value array
cb(sig_array, array(nfunctions*nfunctions))
for j in range(len(function_family)*len(function_family)):
cb(sig_array[j], 0)
for i in range(len(function_family)):
cb(sig_array[i*(nfunctions+1)], sigma[i]**-1)
cb(intermed, matmul(v, sig_array, nfunctions, nfunctions))
cb(ainv, matmul(intermed, ut, nfunctions, nfunctions))
cb(a_coeffs, matmul(ainv, coeff_rhs, nfunctions, 1))
return _linear_comb(
[a_coeffs[ii] for ii in range(hist_len)],
hist_vars)
# }}}
else:
# {{{ static time step
hist_len = len(hist_vars)
nfunctions = len(function_family)
vdm_t = np.zeros((nfunctions, hist_len))
coeff_rhs = np.zeros(nfunctions)
for i in range(nfunctions):
for j in range(hist_len):
vdm_t[i, j] = function_family.evaluate(i, time_values[j])
coeff_rhs[i] = rhs_func(i)
if hist_len == nfunctions:
a_coeffs = la.solve(vdm_t, coeff_rhs)
else:
# SVD-based least squares solve
u, sigma, v = la.svd(vdm_t, full_matrices=False)
ainv = np.dot(v.transpose(), np.dot(la.inv(np.diag(sigma)),
u.transpose()))
a_coeffs = np.dot(ainv, coeff_rhs)
return _linear_comb(a_coeffs, hist_vars)
# }}}
def emit_adams_integration(cb, name_gen,
function_family, time_values, hist_vars, t_start, t_end):
return _emit_func_family_operation(
cb, name_gen, function_family, time_values, hist_vars,
lambda i: (
function_family.antiderivative(i, t_end)
- function_family.antiderivative(i, t_start)))
def emit_adams_extrapolation(cb, name_gen,
function_family, time_values, hist_vars, t_eval):
return _emit_func_family_operation(
cb, name_gen, function_family, time_values, hist_vars,
lambda i: function_family.evaluate(i, t_eval))
# }}}
# {{{ ab method
class AdamsBashforthMethodBuilder(MethodBuilder):
"""
User-supplied context:
<state> + component_id: The value that is integrated
<func> + component_id: The right hand side
.. automethod:: __init__
.. automethod:: generate
"""
def __init__(self, component_id, function_family=None, state_filter_name=None,
hist_length=None, static_dt=False, order=None):
"""
:arg function_family: Accepts an instance of
:class:`AdamsIntegrationFunctionFamily`
or an integer, in which case the classical monomial function family
with the order given by the integer is used.
:arg static_dt: If *True*, changing the timestep during time integration
is not allowed.
"""
if function_family is not None and order is not None:
raise ValueError("may not specify both function_family and order")
if function_family is None:
function_family = order
del order
if isinstance(function_family, int):
function_family = AdamsMonomialIntegrationFunctionFamily(function_family)
super().__init__()
self.function_family = function_family
if hist_length is None:
hist_length = len(function_family)
self.hist_length = hist_length
self.static_dt = static_dt
self.component_id = component_id
# Declare variables
self.step = var("<p>step")
self.function = var("<func>" + component_id)
self.history = \
[var("<p>f_n_minus_" + str(i)) for i in range(hist_length - 1, 0, -1)]
if not self.static_dt:
self.time_history = [
var("<p>t_n_minus_" + str(i))
for i in range(hist_length - 1, 0, -1)]
self.state = var("<state>" + component_id)
self.t = var("<t>")
self.dt = var("<dt>")
if state_filter_name is not None:
self.state_filter = var("<func>" + state_filter_name)
else:
self.state_filter = None
def generate(self):
"""
:returns: :class:`dagrt.language.DAGCode`
"""
from pytools import UniqueNameGenerator
name_gen = UniqueNameGenerator()
from dagrt.language import CodeBuilder, DAGCode
array = var("<builtin>array")
rhs_var = var("rhs_var")
# Initialization
with CodeBuilder(name="initialization") as cb_init:
cb_init(self.step, 1)
# Primary
with CodeBuilder(name="primary") as cb_primary:
if not self.static_dt:
time_history_data = self.time_history + [self.t]
time_hist_var = var(name_gen("time_history"))
cb_primary(time_hist_var, array(self.hist_length))
for i in range(self.hist_length):
cb_primary(time_hist_var[i], time_history_data[i] - self.t)
time_hist = time_hist_var
t_end = self.dt
dt_factor = 1
else:
time_hist = list(range(-self.hist_length+1, 0+1)) # noqa pylint:disable=invalid-unary-operand-type
dt_factor = self.dt
t_end = 1
cb_primary(rhs_var, self.eval_rhs(self.t, self.state))
history = self.history + [rhs_var]
ab_sum = emit_adams_integration(
cb_primary, name_gen,
self.function_family,
time_hist, history,
0, t_end)
state_est = self.state + dt_factor * ab_sum
if self.state_filter is not None:
state_est = self.state_filter(state_est)
cb_primary(self.state, state_est)
# Rotate history and time history.
for i in range(self.hist_length - 1):
cb_primary(self.history[i], history[i + 1])
if not self.static_dt:
cb_primary(self.time_history[i], time_history_data[i + 1])
cb_primary(self.t, self.t + self.dt)
cb_primary.yield_state(expression=self.state,
component_id=self.component_id,
time_id="", time=self.t)
if self.hist_length == 1:
# The first order method requires no bootstrapping.
return DAGCode(
phases={
"initial": cb_init.as_execution_phase(next_phase="primary"),
"primary": cb_primary.as_execution_phase(next_phase="primary")
},
initial_phase="initial")
# Bootstrap
with CodeBuilder(name="bootstrap") as cb_bootstrap:
self.rk_bootstrap(cb_bootstrap)
cb_bootstrap(self.t, self.t + self.dt)
cb_bootstrap.yield_state(expression=self.state,
component_id=self.component_id,
time_id="", time=self.t)
cb_bootstrap(self.step, self.step + 1)
with cb_bootstrap.if_(self.step, "==", self.hist_length):
cb_bootstrap.switch_phase("primary")
return DAGCode(
phases={
"initialization": cb_init.as_execution_phase("bootstrap"),
"bootstrap": cb_bootstrap.as_execution_phase("bootstrap"),
"primary": cb_primary.as_execution_phase("primary"),
},
initial_phase="initialization")
def eval_rhs(self, t, y):
"""Return a node that evaluates the RHS at the given time and
component value."""
from pymbolic.primitives import CallWithKwargs
return CallWithKwargs(function=self.function,
parameters=(),
kw_parameters={"t": t, self.component_id: y})
def rk_bootstrap(self, cb):
"""Initialize the timestepper with an RK method."""
rhs_var = var("rhs_var")
cb(rhs_var, self.eval_rhs(self.t, self.state))
# Save the current RHS to the AB history
for i in range(len(self.history)):
with cb.if_(self.step, "==", i + 1):
cb(self.history[i], rhs_var)
if not self.static_dt:
cb(self.time_history[i], self.t)
from leap.rk import ORDER_TO_RK_METHOD_BUILDER
rk_method = ORDER_TO_RK_METHOD_BUILDER[self.function_family.order]
rk_tableau = tuple(zip(rk_method.c, rk_method.a_explicit))
rk_coeffs = rk_method.output_coeffs
# Stage loop (taken from EmbeddedButcherTableauMethodBuilder)
rhss = [var("rk_rhs_" + str(i)) for i in range(len(rk_tableau))]
for stage_num, (c, coeffs) in enumerate(rk_tableau):
if len(coeffs) == 0:
assert c == 0
cb(rhss[stage_num], rhs_var)
else:
stage = self.state + sum(self.dt * coeff * rhss[j]
for (j, coeff)
in enumerate(coeffs))
if self.state_filter is not None:
stage = self.state_filter(stage)
cb(rhss[stage_num], self.eval_rhs(self.t + c * self.dt, stage))
# Merge the values of the RHSs.
rk_comb = sum(coeff * rhss[j] for j, coeff in enumerate(rk_coeffs))
state_est = self.state + self.dt * rk_comb
if self.state_filter is not None:
state_est = self.state_filter(state_est)
# Assign the value of the new state.
cb(self.state, state_est)
# }}}
# vim: fdm=marker