-
Notifications
You must be signed in to change notification settings - Fork 1.2k
/
implicit_euler_integrator.h
172 lines (147 loc) · 7.09 KB
/
implicit_euler_integrator.h
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
#pragma once
#include <memory>
#include "drake/common/default_scalars.h"
#include "drake/common/drake_copyable.h"
#include "drake/systems/analysis/implicit_integrator.h"
#include "drake/systems/analysis/runge_kutta2_integrator.h"
namespace drake {
namespace systems {
/**
* A first-order, fully implicit integrator with second order error estimation.
*
* This integrator uses the following update rule:<pre>
* x(t+h) = x(t) + h f(t+h,x(t+h))
* </pre>
* where x are the state variables, h is the integration step size, and
* f() returns the time derivatives of the state variables. Contrast this
* update rule to that of an explicit first-order integrator:<pre>
* x(t+h) = x(t) + h f(t, x(t))
* </pre>
* Thus implicit first-order integration must solve a nonlinear system of
* equations to determine *both* the state at t+h and the time derivatives
* of that state at that time. Cast as a nonlinear system of equations,
* we seek the solution to:<pre>
* x(t+h) - x(t) - h f(t+h,x(t+h)) = 0
* </pre>
* given unknowns x(t+h).
*
* This "implicit Euler" method is known to be L-Stable, meaning both that
* applying it at a fixed integration step to the "test" equation `y(t) = eᵏᵗ`
* yields zero (for `k < 0` and `t → ∞`) *and* that it is also A-Stable.
* A-Stability, in turn, means that the method can integrate the linear constant
* coefficient system `dx/dt = Ax` at any step size without the solution
* becoming unstable (growing without bound). The practical effect of
* L-Stability is that the integrator tends to be stable for any given step size
* on an arbitrary system of ordinary differential equations. See
* [Lambert, 1991], Ch. 6 for an approachable discussion on stiff differential
* equations and L- and A-Stability.
*
* This implementation uses Newton-Raphson (NR) and relies upon the obvious
* convergence to a solution for `g = 0` where
* `g(x(t+h)) ≡ x(t+h) - x(t) - h f(t+h,x(t+h))` as `h` becomes sufficiently
* small. It also uses the implicit trapezoid method- fed the result from
* implicit Euler for (hopefully) faster convergence- to compute the error
* estimate. General implementational details were gleaned from [Hairer, 1996].
*
* - [Hairer, 1996] E. Hairer and G. Wanner. Solving Ordinary Differential
* Equations II (Stiff and Differential-Algebraic Problems).
* Springer, 1996.
* - [Lambert, 1991] J. D. Lambert. Numerical Methods for Ordinary Differential
* Equations. John Wiley & Sons, 1991.
*
* @note This integrator uses the integrator accuracy setting, even when run
* in fixed-step mode, to limit the error in the underlying Newton-Raphson
* process. See IntegratorBase::set_target_accuracy() for more info.
* @see ImplicitIntegrator class documentation for information about implicit
* integration methods in general.
*
* @tparam_nonsymbolic_scalar
*/
template <class T>
class ImplicitEulerIntegrator final : public ImplicitIntegrator<T> {
public:
DRAKE_NO_COPY_NO_MOVE_NO_ASSIGN(ImplicitEulerIntegrator)
~ImplicitEulerIntegrator() override = default;
explicit ImplicitEulerIntegrator(const System<T>& system,
Context<T>* context = nullptr)
: ImplicitIntegrator<T>(system, context) {}
/// The integrator supports error estimation.
bool supports_error_estimation() const final { return true; }
/// This first-order integrator uses embedded second order methods to compute
/// estimates of the local truncation error. The order of the asymptotic
/// difference between these two methods (from which the error estimate is
/// computed) is O(h²).
int get_error_estimate_order() const final { return 2; }
private:
int64_t do_get_num_newton_raphson_iterations() const final {
return num_nr_iterations_;
}
int64_t do_get_num_error_estimator_derivative_evaluations() const final {
return num_err_est_function_evaluations_;
}
int64_t do_get_num_error_estimator_derivative_evaluations_for_jacobian()
const final {
return num_err_est_jacobian_function_evaluations_;
}
int64_t do_get_num_error_estimator_newton_raphson_iterations()
const final {
return num_err_est_nr_iterations_;
}
int64_t do_get_num_error_estimator_jacobian_evaluations() const final {
return num_err_est_jacobian_reforms_;
}
int64_t do_get_num_error_estimator_iteration_matrix_factorizations()
const final {
return num_err_est_iter_factorizations_;
}
void DoResetImplicitIntegratorStatistics() final;
static void ComputeAndFactorImplicitEulerIterationMatrix(
const MatrixX<T>& J, const T& h,
typename ImplicitIntegrator<T>::IterationMatrix* iteration_matrix);
static void ComputeAndFactorImplicitTrapezoidIterationMatrix(
const MatrixX<T>& J, const T& h,
typename ImplicitIntegrator<T>::IterationMatrix* iteration_matrix);
void DoInitialize() final;
bool AttemptStepPaired(const T& t0, const T& h, const VectorX<T>& xt0,
VectorX<T>* xtplus_euler, VectorX<T>* xtplus_trap);
bool StepAbstract(const T& t0, const T& h, const VectorX<T>& xt0,
const std::function<VectorX<T>()>& g,
const std::function<
void(const MatrixX<T>&, const T&,
typename ImplicitIntegrator<T>::IterationMatrix*)>&
compute_and_factor_iteration_matrix,
const VectorX<T>& xtplus_guess,
typename ImplicitIntegrator<T>::IterationMatrix*,
VectorX<T>* xtplus, int trial = 1);
bool DoImplicitIntegratorStep(const T& h) final;
bool StepImplicitEuler(const T& t0, const T& h, const VectorX<T>& xt0,
VectorX<T>* xtplus);
bool StepImplicitTrapezoid(const T& t0, const T& h, const VectorX<T>& xt0,
const VectorX<T>& dx0, const VectorX<T>& xtplus_ie, VectorX<T>* xtplus);
// The last computed iteration matrix and factorization for implicit Euler.
typename ImplicitIntegrator<T>::IterationMatrix ie_iteration_matrix_;
// The last computed iteration matrix and factorization for implicit
// trapezoid.
typename ImplicitIntegrator<T>::IterationMatrix itr_iteration_matrix_;
// Vector used in error estimate calculations.
VectorX<T> err_est_vec_;
// The continuous state update vector used during Newton-Raphson.
std::unique_ptr<ContinuousState<T>> dx_state_;
// Variables to avoid heap allocations.
VectorX<T> xt0_, xdot_, xtplus_ie_, xtplus_tr_;
// Various statistics.
int64_t num_nr_iterations_{0};
// Second order Runge-Kutta method for estimating the integration error when
// the requested step size lies below the working step size.
std::unique_ptr<RungeKutta2Integrator<T>> rk2_;
// Implicit trapezoid specific statistics.
int64_t num_err_est_jacobian_reforms_{0};
int64_t num_err_est_iter_factorizations_{0};
int64_t num_err_est_function_evaluations_{0};
int64_t num_err_est_jacobian_function_evaluations_{0};
int64_t num_err_est_nr_iterations_{0};
};
} // namespace systems
} // namespace drake
DRAKE_DECLARE_CLASS_TEMPLATE_INSTANTIATIONS_ON_DEFAULT_NONSYMBOLIC_SCALARS(
class drake::systems::ImplicitEulerIntegrator)