-
Notifications
You must be signed in to change notification settings - Fork 3
/
erk.hpp
515 lines (384 loc) · 13.4 KB
/
erk.hpp
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
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
/*
Rehuel: a simple C++ library for solving ODEs
Copyright 2017-2019, Stefan Paquay (stefanpaquay@gmail.com)
Licensed under the Apache License, Version 2.0 (the "License");
you may not use this file except in compliance with the License.
You may obtain a copy of the License at
http://www.apache.org/licenses/LICENSE-2.0
Unless required by applicable law or agreed to in writing, software
distributed under the License is distributed on an "AS IS" BASIS,
WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
See the License for the specific language governing permissions and
limitations under the License.
============================================================================= */
/**
\file erk.hpp
\brief Functions related to performing time integration with
explicit Runge-Kutta (RK) methods.
*/
#ifndef ERK_HPP
#define ERK_HPP
#include <cassert>
#include <limits>
#include <iomanip>
#include "arma_include.hpp"
#include "enums.hpp"
#include "my_timer.hpp"
#include "newton.hpp"
#include "options.hpp"
#include "output.hpp"
/**
\brief Contains functions related to implicit Runge-Kutta methods.
*/
namespace erk {
#ifdef DEBUG_OUTPUT
constexpr const bool debug = true;
#else
constexpr const bool debug = false;
#endif // DEBUG_OUTPUT
typedef arma::mat mat_type;
typedef arma::vec vec_type;
/**
Contains the Butcher tableau plus time step size.
*/
struct solver_coeffs
{
const char *name; ///< Human-friendly name for the method.
vec_type b; ///< weights for the new y-value
vec_type c; ///< these set the intermediate time points
mat_type A; ///< alpha coefficients in Butcher tableau
vec_type b2; ///< weights for the new y-value of the embedded RK method
int order; ///< Local convergence order for main method
int order2; ///< Local convergence order for embedded method
/// If the method satisfies first-same-as-last (FSAL), set this to
/// true to enable the optimizations associated with FSAL.
bool FSAL;
/// This gamma is used as weight for the current value for embedding.
double gamma;
/// This matrix defines the interpolating polynomial, if available.
mat_type b_interp;
};
/**
\brief options for the time integrator.
*/
struct solver_options : common_solver_options {
/// \brief Constructor with default values.
solver_options() : adaptive_step_size(true),
extrapolate_stage(false)
{ }
~solver_options()
{ }
/// If true, attempt to perform adaptive time stepping using
/// an embedded pair.
bool adaptive_step_size;
/// If true, use current stages and extrapolate to the next time level.
/// \note This is hard for RK-methods not based on quadrature, so is
/// generally not supported.
bool extrapolate_stage;
};
static std::map<int,std::string> rk_method_to_string = {
FOREACH_ERK_METHOD(GENERATE_STRING)
};
static std::map<std::string,int> rk_string_to_method = {
FOREACH_ERK_METHOD(GENERATE_MAP)
};
/**
\brief a struct that contains time stamps and stages that can be used for
constructing the solution all time points in the interval (dense output).
*/
struct rk_output : basic_output
{
struct counters {
counters() : attempt(0), reject_err(0), fun_evals(0) {}
std::size_t attempt, reject_err;
std::size_t fun_evals;
};
std::vector<vec_type> stages;
std::vector<vec_type> err_est;
std::vector<double> err;
double elapsed_time, accept_frac;
counters count;
};
/**
\brief Returns a vector with all method names.
*/
std::vector<std::string> all_method_names();
/**
\brief Returns coefficients belonging to the given method.
See erk::rk_methods for all methods.
\note If the method is not recognized, the coefficients
returned will not pass verify_solver_coefficients.
\param method The method to return coefficients for.
\returns coefficients belonging to given method.
*/
solver_coeffs get_coefficients( int method );
/**
\brief Default solver options:
*/
solver_options default_solver_options();
/**
\brief Checks if all options are set to sane values.
\returns true if all options checked out, false otherwise.
*/
bool verify_solver_options( const solver_options &opts );
/**
\brief Checks whether or not the given coefficients are consistent in size.
\param sc the coefficients to check.
\returns true if the coefficients are valid, false otherwise.
*/
bool verify_solver_coeffs( const erk::solver_coeffs &sc );
/**
\brief Converts a string with a method name to an int.
\param name A string describing the method.
\returns the enum corresponding to given method. See \ref rk_methods
*/
int name_to_method( const std::string &name );
/**
\brief Converts method code to a human-readable string.
\note This output shall satisfy
name_to_method( method_to_name( method ) ) == method.
\param method The method to convert to a name.
\returns a string literal representing the method.
*/
const char *method_to_name( int method );
// This is copy-pasta from irk.cpp, should be moved I think.
vec_type project_b( double theta, const erk::solver_coeffs &sc );
/**
\brief Takes stage matrix and assigns the last to the first.
\param Ks stage matrix
\param Ns number of stages.
*/
void apply_fsal(mat_type &Ks, std::size_t Ns)
{
Ks.col(0) = std::move(Ks.col(Ns-1));
}
/**
\brief Takes stage matrix does nothing.
\param Ks stage matrix
\param Ns number of stages.
*/
void no_apply_fsal_dummy(mat_type &Ks, std::size_t Ns)
{ }
/**
\brief Guts of the explicit RK integrator.
Time-integrates a given ODE from t0 to t1, starting at y0
t_vals and y_vals shall be unmodified upon failure.
\param func Functor of the ODE to integrate
\param t0 Starting time
\param t1 Final time
\param y0 Initial values
\param solver_opts Options for the internal solver.
\param dt Initial time step size.
\param sc Coefficients of the solver.
\returns an output struct with the solution.
*/
template <typename functor_type> inline
rk_output erk_guts(functor_type &func, double t0, double t1, const vec_type &y0,
const solver_options &solver_opts, double dt,
const solver_coeffs &sc )
{
if( t0 + dt > t1 ){
std::cerr << " Rehuel: Initial dt (" << dt;
dt = t1 - t0;
std::cerr << ") too large for interval! Reducing to "
<< dt << "\n";
}
std::cerr << " Rehuel: Integrating over interval [ "
<< t0 << ", " << t1 << " ]...\n"
<< " Method = " << sc.name << "\n";
// Explicit RK methods are a lot simpler.
// First you compute the k stages explicitly:
my_timer timer;
double t = t0;
rk_output sol;
sol.status = SUCCESS;
assert(dt > 0 && "Cannot use time step size <= 0!");
std::size_t Neq = y0.size();
std::size_t Ns = sc.b.size();
vec_type y = y0;
// Ks is the stages at the new time step.
mat_type Ks(Neq,Ns);
long long int step = 0;
// For time step size control.
double dts[3] = {dt, dt, dt}, errs[3] = {0.9,0.9,0.9};
if( solver_opts.out_interval > 0 ){
std::cerr << " Rehuel: step t dt err\n";
}
double err = 0.0;
vec_type err_est = arma::zeros(Neq);
sol.t_vals.push_back(t);
sol.y_vals.push_back(y);
sol.stages.push_back(vectorise(Ks));
sol.err_est.push_back( err_est );
sol.err.push_back( err );
// This will keep track of error estimates during integration.
std::size_t min_order = std::min( sc.order, sc.order2 );
// If your method has FSAL, you never have to compute the first stage
// after the first step.
std::size_t stage_iter_start = 0;
auto fsal_hook_fptr = no_apply_fsal_dummy;
// By wrapping the function call in this lambda, you can more easily
// count the number of function evaluations.
auto eval_fun = [&func,&sol](double t, const vec_type &Y)
{ ++sol.count.fun_evals; return func.fun(t, Y); };
if (sc.FSAL) {
stage_iter_start = 1;
Ks.col(0) = eval_fun(t, y0);
fsal_hook_fptr = apply_fsal;
}
while( t < t1 ) {
// **************** Calculate stages: ************
// Make sure you stop exactly at t = t1.
if( t + dt > t1 ){
dt = t1 - t;
}
sol.count.attempt++;
if (solver_opts.max_steps >= 0 &&
step > solver_opts.max_steps) {
std::cerr << " Rehuel: Maximum number of attempts exceeded.\n";
sol.status = ERROR_MAX_STEPS_EXCEEDED;
return sol;
}
int integrator_status = 0;
// Formula for explicit stages are
// k_i = f(t + ci*dt, y0 + sum_{j=1}^{i-1} A(i,j)*k_j)
for (std::size_t i = stage_iter_start; i < Ns; ++i) {
vec_type tmp = y;
for (std::size_t j = 0; j < i; ++j) {
tmp += dt*sc.A(i,j)*Ks.col(j);
}
Ks.col(i) = eval_fun(t + sc.c(i)*dt, tmp);
}
// ************* Form solution at t + dt: ***********
vec_type delta_y = Ks*sc.b;
vec_type y_n = y + dt*delta_y;
double new_dt = dt;
// If you have no adaptive step size, error calculation
// might not be very sensible.
if (solver_opts.adaptive_step_size) {
vec_type delta_alt = Ks*sc.b2;
// ************* Error estimate: ***********
double err_tot = 0.0;
double atol = solver_opts.abs_tol, rtol = solver_opts.rel_tol;
// err_est is ||y1 - yhat1|| in Wanner & Hairer.
err_est = dt*(delta_alt - delta_y);
double n = Neq;
for (std::size_t i = 0; i < Neq; ++i) {
double erri = err_est(i);
double y0i = std::fabs(y(i));
double y1i = std::fabs(y_n(i));
double sci = atol + rtol * std::max(y0i, y1i);
double add = erri / sci;
err_tot += add*add;
}
err = std::sqrt(err_tot / n);
assert( err_tot >= 0.0 && "Error cannot be negative!" );
if (err < machine_precision) {
err = machine_precision;
}
errs[2] = errs[1];
errs[1] = errs[0];
errs[0] = err;
// ************* Adaptive time step size control: ***********
// TODO: This might need optimization geared
// to explicit methods.
// Error is too large to tolerate:
if (solver_opts.adaptive_step_size && (err >= 1.0)) {
integrator_status = 1;
sol.count.reject_err++;
}
double fac = 0.9;
double expt = 1.0 / (1.0 + min_order);
double err_inv = 1.0 / err;
double scale_27 = std::pow(err_inv, expt);
// double dt_rat = dts[0] / dts[1];
// double err_frac = errs[1] / errs[0];
// if (errs[1] == 0 || errs[0] == 0) {
// err_frac = 1.0;
// }
// double err_rat = std::pow(err_frac, expt);
double min_scale = scale_27; // std::min(scale_27, scale_28);
new_dt = fac * dt * std::min(4.0, min_scale);
if (solver_opts.max_dt > 0) {
new_dt = std::min(solver_opts.max_dt, new_dt);
}
}
// ********************* Output if user requested **************
if (solver_opts.out_interval > 0 &&
(step % solver_opts.out_interval == 0) ) {
std::cerr << " Rehuel: " << step << " " << t
<< " " << dt << " " << err << "\n";
}
// ********************* Update y and time ***************
if (!solver_opts.adaptive_step_size || integrator_status == 0) {
y = y_n;
t += dt;
++step;
sol.t_vals.push_back(t);
sol.y_vals.push_back(y_n);
// Since K is a matrix, it needs to be flattened:
sol.stages.push_back(arma::vectorise(Ks));
sol.err_est.push_back(err_est);
sol.err.push_back(err);
// If you reach here, your new time step has been
// accepted and your last stage can now be uesd
// as your first stage:
fsal_hook_fptr(Ks, Ns);
}
// **************** Set the new time step size. *********************
if (solver_opts.adaptive_step_size) {
dt = new_dt;
}
dts[2] = dts[1];
dts[1] = dts[0];
dts[0] = dt;
}
double elapsed = timer.toc();
sol.elapsed_time = elapsed;
sol.accept_frac = static_cast<double>(step) / sol.count.attempt;
return sol;
}
/**
\brief Time-integrate a given ODE from t0 to t1, starting at y0
t_vals and y_vals shall be unmodified upon failure.
\param func Functor of the ODE to integrate
\param t0 Starting time
\param t1 Final time
\param y0 Initial values
\param dt Initial time step size.
\param solver_opts Options for the internal solver.
\returns a status code (see \ref odeint_status_codes)
*/
template <typename functor_type> inline
rk_output odeint(functor_type &func, double t0, double t1, const vec_type &y0,
solver_options solver_opts,
int method = erk::DORMAND_PRINCE_54, double dt = 1e-6)
{
solver_coeffs sc = get_coefficients(method);
if (solver_opts.adaptive_step_size && sc.b2.size() == 0) {
std::cerr << " Rehuel: WARNING: Cannot have adaptive time "
<< "step with non-embedding method! Disabling "
<< "adaptive time step size!\n";
solver_opts.adaptive_step_size = false;
}
assert( verify_solver_coeffs( sc ) && "Invalid solver coefficients!" );
return erk_guts(func, t0, t1, y0, solver_opts, dt, sc);
}
/**
\brief Time-integrate a given ODE from t0 to t1, starting at y0.
This function is supposed to provide a sane "default" explicit solver,
a la ode45 in Matlab.
\param func Functor of the ODE to integrate
\param t0 Starting time
\param t1 Final time
\param y0 Initial values
\returns a struct with the solution and info about the solution quality.
*/
template <typename functor_type> inline
rk_output odeint( functor_type &func, double t0, double t1, const vec_type &y0)
{
solver_options s_opts = default_solver_options();
return odeint(func, t0, t1, y0, s_opts);
}
} // namespace erk
#endif // ERK_HPP