13#ifndef SANDALS_RUNGEKUTTA_HXX
14#define SANDALS_RUNGEKUTTA_HXX
37 template <Integer S, Integer N, Integer M>
41 using MatrixK = Eigen::Matrix<Real, N, S>;
42 using MatrixJ = Eigen::Matrix<Real, N*S, N*S>;
43 using VectorP = Eigen::Matrix<Real, N+M, 1>;
44 using MatrixP = Eigen::Matrix<Real, N+M, N+M>;
45 using NewtonX = Optimist::RootFinder::Newton<N>;
46 using NewtonK = Optimist::RootFinder::Newton<N*S>;
57 using Time = Eigen::Vector<Real, Eigen::Dynamic>;
75 Eigen::FullPivLU<MatrixP>
m_lu;
121 bool is_erk()
const {
return this->m_tableau.
type == Type::ERK;}
127 bool is_irk()
const {
return this->m_tableau.
type == Type::IRK;}
157 std::string
name()
const {
return this->m_tableau.
name;}
331 this->m_newtonX.verbose_mode(t_verbose);
332 this->m_newtonK.verbose_mode(t_verbose);
391 {this->m_max_projection_iterationsations = t_max_projection_iterations;}
452 std::max(x.array().abs().maxCoeff(), x_e.array().abs().maxCoeff())};
453 Real truncation_error{(x - x_e).array().abs().maxCoeff()};
456 1.0/std::max(this->m_tableau.
order, this->m_tableau.order_e))));
464 std::ostringstream os;
466 <<
"Runge-Kutta method:\t" << this->
name() << std::endl
467 <<
"\t- order:\t" << this->
order() << std::endl
468 <<
"\t- stages:\t" << this->
stages() << std::endl
470 switch (this->
type()) {
471 case Type::ERK: os <<
"explicit";
break;
472 case Type::IRK: os <<
"implicit";
break;
473 case Type::DIRK: os <<
"diagonally implicit";
break;
477 <<
"\t- embedded:\t" << this->
is_embedded() << std::endl;
479 os <<
"\t- system:\t" << this->m_system->name() << std::endl;
481 os <<
"\t- system:\t" <<
"none" << std::endl;
522 for (
Integer i{0}; i < S; ++i) {
523 x_node = x_old + K(all, seqN(0, i)) * this->m_tableau.
A(i, seqN(0, i)).transpose();
525 K.col(i) = h_old *
static_cast<Explicit<N, M> const *
>(this->m_system.get())->
f(x_node, t_old + h_old*this->m_tableau.
c(i));
527 K.col(i) = h_old *
static_cast<Explicit<N, M> const *
>(this->m_system.get())->
f_reverse(x_node, t_old + h_old*this->m_tableau.
c(i));
530 if (!K.allFinite()) {
return false;}
533 x_new = x_old + K * this->m_tableau.
b;
537 VectorN x_emb = x_old + K * this->m_tableau.
b_e;
567 VectorN x_node(x + K(all, seqN(0, s)) * this->m_tableau.
A(s, seqN(0, s)).transpose());
569 fun = this->m_system->F(x_node, K.col(s)/h, t + h * this->m_tableau.c(s));
571 fun = this->m_system->F_reverse(x_node, K.col(s)/h, t + h * this->m_tableau.c(s));
610 VectorN x_node(x + K(all, seqN(0, s)) * this->m_tableau.
A(s, seqN(0, s)).transpose());
612 jac = this->m_system->JF_x_dot(x_node, K.col(s)/h, t + h * this->m_tableau.c(s)) / h;
614 jac = this->m_system->JF_x_dot_reverse(x_node, K.col(s)/h, t + h * this->m_tableau.c(s)) / h;
635 VectorN K_ini(VectorN::Zero());
638 for (
Integer s{0}; s < S; ++s) {
639 if (this->m_newtonX.solve(
640 [
this, s, &K, &x_old, t_old, h_old](
VectorN const &K_fun,
VectorN &fun)
641 {K.col(s) = K_fun; this->erk_implicit_function(s, x_old, t_old, h_old, K, fun);},
642 [
this, s, &K, &x_old, t_old, h_old](
VectorN const &K_jac,
MatrixN &jac)
643 {K.col(s) = K_jac; this->erk_implicit_jacobian(s, x_old, t_old, h_old, K, jac);},
652 x_new = x_old + K * this->m_tableau.
b;
656 VectorN x_emb(x_old + K * this->m_tableau.
b_e);
698 MatrixK K_mat{K.reshaped(N, S)};
700 for (
Integer i{0}; i < S; ++i) {
701 x_node = x + K_mat * this->m_tableau.
A.row(i).transpose();
703 fun_mat.col(i) = this->m_system->F(x_node, K_mat.col(i)/h, t + h * this->m_tableau.c(i));
705 fun_mat.col(i) = this->m_system->F_reverse(x_node, K_mat.col(i)/h, t + h * this->m_tableau.c(i));
708 fun = fun_mat.reshaped(N*S, 1);
755 MatrixK K_mat{K.reshaped(N, S)};
759 auto idx = seqN(0, N), jdx = seqN(0, N);
760 for (
Integer i{0}; i < S; ++i) {
761 t_node = t + h * this->m_tableau.
c(i);
762 x_node = x + K_mat * this->m_tableau.
A.row(i).transpose();
765 x_dot_node = K_mat.col(i) / h;
767 JF_x = this->m_system->JF_x(x_node, x_dot_node, t_node);
768 JF_x_dot = this->m_system->JF_x_dot(x_node, x_dot_node, t_node);
770 JF_x = this->m_system->JF_x_reverse(x_node, x_dot_node, t_node);
771 JF_x_dot = this->m_system->JF_x_dot_reverse(x_node, x_dot_node, t_node);
776 for (
Integer j{0}; j < S; ++j) {
779 jac(idx, jdx) = this->m_tableau.
A(i,j) * JF_x + JF_x_dot / h;
781 jac(idx, jdx) = this->m_tableau.
A(i,j) * JF_x;
803 VectorK K_ini(VectorK::Zero());
806 if (!this->m_newtonK.solve(
808 {this->irk_function(x_old, t_old, h_old, K_fun, fun);},
810 {this->irk_jacobian(x_old, t_old, h_old, K_jac, jac);},
815 x_new = x_old + K.reshaped(N, S) * this->m_tableau.
b;
819 VectorN x_emb(x_old + K.reshaped(N, S) * this->m_tableau.b_e);
858 VectorN x_node(x + K(all, seqN(0, n+1)) * this->m_tableau.
A(n, seqN(0, n+1)).transpose());
860 fun = this->m_system->F(x_node, K.col(n)/h, t + h * this->m_tableau.c(n));
862 fun = this->m_system->F_reverse(x_node, K.col(n)/h, t + h * this->m_tableau.c(n));
903 Real t_node{t + h * this->m_tableau.
c(n)};
904 VectorN x_node(x + K(all, seqN(0, n+1)) * this->m_tableau.
A(n, seqN(0, n+1)).transpose());
905 VectorN x_dot_node(K.col(n)/h);
907 jac = this->m_tableau.
A(n,n) * this->m_system->JF_x(x_node, x_dot_node, t_node) +
908 this->m_system->JF_x_dot(x_node, x_dot_node, t_node) / h;
910 jac = this->m_tableau.
A(n,n) * this->m_system->JF_x_reverse(x_node, x_dot_node, t_node) +
911 this->m_system->JF_x_dot_reverse(x_node, x_dot_node, t_node) / h;
932 VectorN K_ini(VectorN::Zero());
935 for (
Integer n{0}; n < S; ++n) {
936 if (this->m_newtonX.solve(
937 [
this, n, &K, &x_old, t_old, h_old](
VectorN const &K_fun,
VectorN &fun)
938 {K.col(n) = K_fun; this->dirk_function(n, x_old, t_old, h_old, K, fun);},
939 [
this, n, &K, &x_old, t_old, h_old](
VectorN const &K_jac,
MatrixN &jac)
940 {K.col(n) = K_jac; this->dirk_jacobian(n, x_old, t_old, h_old, K, jac);},
949 x_new = x_old + K * this->m_tableau.
b;
953 VectorN x_emb(x_old + K * this->m_tableau.
b_e);
972 #define CMD "Sandals::RungeKutta::step(...): "
974 SANDALS_ASSERT(this->m_system->in_domain(x_old, t_old),
CMD "in " << this->m_tableau.name <<
975 " solver, at t = " << t_old <<
", x = " << x_old.transpose() <<
", system out of domain.");
977 if (this->
is_erk() && this->m_system->is_explicit()) {
979 }
else if (this->
is_erk() && this->m_system->is_implicit()) {
982 return this->
dirk_step(x_old, t_old, h_old, x_new, h_new);
984 return this->
irk_step(x_old, t_old, h_old, x_new, h_new);
1005 #define CMD "Sandals::RungeKutta::advance(...): "
1009 h_old <<
", expected > 0.");
1012 if (!this->
step(x_old, t_old, h_old, x_new, h_new))
1015 Real t_tmp{t_old}, h_tmp{h_old /
Real(2.0)};
1022 if (this->
step(x_tmp, t_tmp, h_tmp, x_new, h_new_tmp)) {
1028 if (k > 0 && k < max_k) {
1032 h_tmp =
Real(2.0) * h_tmp;
1035 ", integration succedded disable one substepping layer.");
1041 SANDALS_ASSERT(std::isfinite(x_tmp.maxCoeff()),
CMD "in " << this->m_tableau.name <<
1042 " solver, at t = " << t_tmp <<
", ||x||_inf = inf, computation interrupted.");
1049 t_tmp <<
", integration failed with h = " << h_tmp <<
", aborting.");
1054 "at t = " << t_tmp <<
", integration failed, adding substepping layer.");}
1071 if (this->
project(x_new, t_old + h_new, x_projected)) {
1072 x_new = x_projected;
1098 sol.
resize(t_mesh.size());
1101 sol.
t(0) = t_mesh(0);
1103 sol.
h.col(0) = this->m_system->h(ics, t_mesh(0));
1108 Real t_step{t_mesh(0)}, h_step{t_mesh(1) - t_mesh(0)}, h_tmp_step{h_step}, h_new_step;
1109 bool mesh_point_bool, saturation_bool;
1113 if (!this->
advance(sol.
x.col(
step), t_step, h_step, x_step, h_new_step)) {
return false;}
1122 h_tmp_step = h_new_step;
1123 h_step = t_mesh(
step+1) - t_step;
1125 h_step = h_new_step;
1132 h_step = h_tmp_step;
1135 sol.
t(
step) = t_step;
1136 sol.
x.col(
step) = x_step;
1137 sol.
h.col(
step) = this->m_system->h(x_step, t_step);
1140 if (std::abs(t_step - t_mesh(last)) <
SQRT_EPSILON) {
break;}
1160 #define CMD "Sandals::RungeKutta::adaptive_solve(...): "
1165 return this->
solve(t_mesh, ics, sol);
1168 return this->
solve(t_mesh, ics, sol);
1172 Real h_step{t_mesh(1) - t_mesh(0)}, h_new_step, scale{100.0};
1173 Real h_min{std::max(this->
m_min_step, h_step/scale)}, h_max{scale*h_step};
1175 Integer safety_length{
static_cast<Integer>(std::ceil(std::abs(t_mesh(last) - t_mesh(0))/(2.0*h_min)))};
1176 sol.
resize(safety_length);
1178 sol.
resize(t_mesh.size());
1182 sol.
t(0) = t_mesh(0);
1184 sol.
h.col(0) = this->m_system->h(ics, t_mesh(0));
1196 h_step = std::max(std::min(h_new_step, h_max), h_min);
1203 sol.
x.col(
step+1) = x_step;
1204 sol.
h.col(
step+1) = this->m_system->h(x_step, sol.
t(
step+1));
1207 if (sol.
t(
step+1) + h_step > t_mesh(last)) {
break;}
1232 #define CMD "Sandals::RungeKutta::project(...): "
1243 A.template block<N, N>(0, 0) = MatrixN::Identity();
1254 h = this->m_system->h(x_projected, t);
1255 Jh_x = this->m_system->Jh_x(x_projected, t);
1258 if (h.norm() < this->m_projection_tolerance) {
return true;}
1261 A.template block<N, M>(0, N) = Jh_x.transpose();
1262 A.template block<M, N>(N, 0) = Jh_x;
1263 b.template head<N>() = x - x_projected;
1264 b.template tail<M>() = -h;
1267 this->m_lu.compute(
A);
1269 x_step = this->m_lu.solve(
b);
1272 if (x_step.norm() < this->m_projection_tolerance * this->m_projection_tolerance) {
return false;}
1275 x_projected.noalias() += x_step(Eigen::seqN(0, N));
1298 std::vector<Integer>
const & projected_invariants,
VectorN &x_projected)
const
1300 #define CMD "Sandals::RungeKutta::project_ics(...): "
1314 A.block(0, 0, X, X) = MatrixX::Identity(X+H, X+H);
1315 Eigen::FullPivLU<MatrixX> lu;
1326 h = this->m_system->h(x_projected, t);
1327 Jh_x = this->m_system->Jh_x(x_projected, t);
1330 h = h(projected_invariants);
1331 Jh_x = Jh_x(projected_invariants, projected_equations);
1334 if (h.norm() < this->m_projection_tolerance) {
return true;}
1337 A.block(0, X, X, H) = Jh_x.transpose();
1338 A.block(X, 0, H, X) = Jh_x;
1339 b.head(X) = x(projected_equations) - x_projected(projected_equations);
1345 x_step = this->m_lu.solve(
b);
1348 if (x_step.norm() < this->m_projection_tolerance * this->m_projection_tolerance) {
return false;}
1351 x_projected(projected_equations).noalias() += x_step;
1373 #define CMD "Sandals::RungeKutta::estimate_order(...): "
1377 for (
Integer i{0}; i < static_cast<Integer>(t_mesh.size()); ++i) {
1381 CMD "expected the same initial time.");
1383 CMD "expected the same final time.");
1386 for (
Integer j{1}; j < static_cast<Integer>(t_mesh[i].size()); ++j) {
1388 CMD "expected a fixed step.");
1395 VectorX h_vec(t_mesh.size()), e_vec(t_mesh.size());
1396 for (
Integer i{0}; i < static_cast<Integer>(t_mesh.size()); ++i) {
1398 "the" << i <<
"-th time mesh.");
1399 sol_ana = sol(sol_num.
t);
1401 CMD "expected the same number of states in analytical solution.");
1403 CMD "expected the same number of steps in analytical solution.");
1404 h_vec(i) = std::abs(sol_num.
t(1) - sol_num.
t(0));
1405 e_vec(i) = (sol_ana - sol_num.
x).array().abs().maxCoeff();
1411 return ((
A.transpose() *
A).ldlt().solve(
A.transpose() *
b))(0);
#define SANDALS_ASSERT(COND, MSG)
Definition Sandals.hh:43
#define SANDALS_WARNING(MSG)
Definition Sandals.hh:52
Class container for the system of explicit ODEs.
Definition Explicit.hxx:38
VectorF f_reverse(VectorF const &x, Real t) const
Definition Explicit.hxx:156
virtual VectorF f(VectorF const &x, Real t) const =0
Eigen::Vector< Real, N > VectorF
Definition Implicit.hxx:43
Eigen::Matrix< Real, N, N > MatrixJF
Definition Implicit.hxx:44
Eigen::Vector< Real, M > VectorH
Definition Implicit.hxx:45
Eigen::Matrix< Real, M, N > MatrixJH
Definition Implicit.hxx:46
std::shared_ptr< Implicit< N, M > > Pointer
Definition Implicit.hxx:42
bool m_reverse
Definition RungeKutta.hxx:73
Real relative_tolerance()
Definition RungeKutta.hxx:229
System system()
Definition RungeKutta.hxx:199
Tableau< S > m_tableau
Definition RungeKutta.hxx:62
Real m_safety_factor
Definition RungeKutta.hxx:66
void disable_projection()
Definition RungeKutta.hxx:413
void enable_adaptive_mode()
Definition RungeKutta.hxx:312
bool step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hxx:970
void erk_implicit_function(Integer s, VectorN const &x, Real t, Real h, MatrixK const &K, VectorN &fun) const
Definition RungeKutta.hxx:563
void max_substeps(Integer t_max_substeps)
Definition RungeKutta.hxx:295
bool m_verbose
Definition RungeKutta.hxx:72
Eigen::FullPivLU< MatrixP > m_lu
Definition RungeKutta.hxx:75
Real projection_tolerance()
Definition RungeKutta.hxx:371
void min_safety_factor(Real t_min_safety_factor)
Definition RungeKutta.hxx:259
Real absolute_tolerance()
Definition RungeKutta.hxx:217
bool erk_implicit_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hxx:631
void system(System t_system)
Definition RungeKutta.hxx:205
Real m_min_step
Definition RungeKutta.hxx:69
Real m_min_safety_factor
Definition RungeKutta.hxx:67
Integer m_max_projection_iterations
Definition RungeKutta.hxx:77
void relative_tolerance(Real t_relative_tolerance)
Definition RungeKutta.hxx:235
bool is_embedded() const
Definition RungeKutta.hxx:169
void safety_factor(Real t_safety_factor)
Definition RungeKutta.hxx:247
Real estimate_order(std::vector< VectorX > const &t_mesh, VectorN const &ics, std::function< MatrixX(VectorX)> &sol)
Definition RungeKutta.hxx:1369
Integer stages() const
Definition RungeKutta.hxx:151
VectorS b_embedded() const
Definition RungeKutta.hxx:187
RungeKutta(Tableau< S > const &t_tableau)
Definition RungeKutta.hxx:96
Integer m_max_substeps
Definition RungeKutta.hxx:70
Real m_absolute_tolerance
Definition RungeKutta.hxx:64
void projection_tolerance(Real t_projection_tolerance)
Definition RungeKutta.hxx:377
void absolute_tolerance(Real t_absolute_tolerance)
Definition RungeKutta.hxx:223
Eigen::Matrix< Real, N+M, N+M > MatrixP
Definition RungeKutta.hxx:44
Real estimate_step(VectorN const &x, VectorN const &x_e, Real h_k) const
Definition RungeKutta.hxx:449
void projection(bool t_projection)
Definition RungeKutta.hxx:403
bool is_dirk() const
Definition RungeKutta.hxx:133
bool reverse_mode()
Definition RungeKutta.hxx:349
Type type() const
Definition RungeKutta.hxx:115
Eigen::Matrix< Real, N+M, 1 > VectorP
Definition RungeKutta.hxx:43
Integer order() const
Definition RungeKutta.hxx:163
bool m_adaptive
Definition RungeKutta.hxx:71
Tableau< S > & tableau()
Definition RungeKutta.hxx:139
bool dirk_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hxx:928
std::string info() const
Definition RungeKutta.hxx:463
bool is_irk() const
Definition RungeKutta.hxx:127
System m_system
Definition RungeKutta.hxx:63
Integer & max_substeps()
Definition RungeKutta.hxx:289
void disable_verbose_mode()
Definition RungeKutta.hxx:343
void dirk_jacobian(Integer n, VectorN const &x, Real t, Real h, MatrixK const &K, MatrixN &jac) const
Definition RungeKutta.hxx:899
void irk_jacobian(VectorN const &x, Real t, Real h, VectorK const &K, MatrixJ &jac) const
Definition RungeKutta.hxx:747
RungeKutta & operator=(RungeKutta const &)=delete
RungeKutta(const RungeKutta &)=delete
VectorS c() const
Definition RungeKutta.hxx:193
bool adaptive_mode()
Definition RungeKutta.hxx:301
void enable_reverse_mode()
Definition RungeKutta.hxx:360
void disable_adaptive_mode()
Definition RungeKutta.hxx:317
Real & safety_factor()
Definition RungeKutta.hxx:241
typename Implicit< N, M >::MatrixJF MatrixN
Definition RungeKutta.hxx:50
Eigen::Matrix< Real, N *S, N *S > MatrixJ
Definition RungeKutta.hxx:42
typename Tableau< S >::Matrix MatrixS
Definition RungeKutta.hxx:48
std::string name() const
Definition RungeKutta.hxx:157
RungeKutta(Tableau< S > const &t_tableau, System t_system)
Definition RungeKutta.hxx:106
typename Implicit< N, M >::Pointer System
Definition RungeKutta.hxx:55
Real & max_safety_factor()
Definition RungeKutta.hxx:265
void enable_verbose_mode()
Definition RungeKutta.hxx:338
void verbose_mode(bool t_verbose)
Definition RungeKutta.hxx:329
Integer & max_projection_iterations()
Definition RungeKutta.hxx:384
bool verbose_mode()
Definition RungeKutta.hxx:323
Tableau< S > const & tableau() const
Definition RungeKutta.hxx:145
typename Tableau< S >::Type Type
Definition RungeKutta.hxx:56
Eigen::Vector< Real, N *S > VectorK
Definition RungeKutta.hxx:40
bool m_projection
Definition RungeKutta.hxx:78
typename Implicit< N, M >::VectorH VectorM
Definition RungeKutta.hxx:51
Eigen::Matrix< Real, N, S > MatrixK
Definition RungeKutta.hxx:41
typename Implicit< N, M >::MatrixJH MatrixM
Definition RungeKutta.hxx:52
bool erk_explicit_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new) const
Definition RungeKutta.hxx:514
typename Implicit< N, M >::VectorF VectorN
Definition RungeKutta.hxx:49
Real m_relative_tolerance
Definition RungeKutta.hxx:65
bool project_ics(VectorN const &x, Real t, std::vector< Integer > const &projected_equations, std::vector< Integer > const &projected_invariants, VectorN &x_projected) const
Definition RungeKutta.hxx:1297
void min_step(Real t_min_step)
Definition RungeKutta.hxx:283
bool is_erk() const
Definition RungeKutta.hxx:121
void reverse(bool t_reverse)
Definition RungeKutta.hxx:355
typename Tableau< S >::Vector VectorS
Definition RungeKutta.hxx:47
VectorS b() const
Definition RungeKutta.hxx:181
void dirk_function(Integer n, VectorN const &x, Real t, Real h, MatrixK const &K, VectorN &fun) const
Definition RungeKutta.hxx:854
void irk_function(VectorN const &x, Real t, Real h, VectorK const &K, VectorK &fun) const
Definition RungeKutta.hxx:695
bool solve(VectorX const &t_mesh, VectorN const &ics, Solution< N, M > &sol)
Definition RungeKutta.hxx:1093
void adaptive(bool t_adaptive)
Definition RungeKutta.hxx:307
bool irk_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hxx:800
bool projection()
Definition RungeKutta.hxx:397
NewtonK m_newtonK
Definition RungeKutta.hxx:61
void enable_projection()
Definition RungeKutta.hxx:408
bool project(VectorN const &x, Real t, VectorN &x_projected)
Definition RungeKutta.hxx:1230
void max_projection_iterations(Integer t_max_projection_iterations)
Definition RungeKutta.hxx:390
Real m_projection_tolerance
Definition RungeKutta.hxx:76
void disable_reverse_mode()
Definition RungeKutta.hxx:365
void max_safety_factor(Real t_max_safety_factor)
Definition RungeKutta.hxx:271
bool has_system()
Definition RungeKutta.hxx:211
bool advance(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hxx:1003
Real m_max_safety_factor
Definition RungeKutta.hxx:68
Eigen::Vector< Real, Eigen::Dynamic > Time
Definition RungeKutta.hxx:57
bool adaptive_solve(VectorX const &t_mesh, VectorN const &ics, Solution< N, M > &sol)
Definition RungeKutta.hxx:1155
Real & min_step()
Definition RungeKutta.hxx:277
Real & min_safety_factor()
Definition RungeKutta.hxx:253
Optimist::RootFinder::Newton< N *S > NewtonK
Definition RungeKutta.hxx:46
void erk_implicit_jacobian(Integer s, VectorN const &x, Real t, Real h, MatrixK const &K, MatrixN &jac) const
Definition RungeKutta.hxx:606
Optimist::RootFinder::Newton< N > NewtonX
Definition RungeKutta.hxx:45
MatrixS A() const
Definition RungeKutta.hxx:175
NewtonX m_newtonX
Definition RungeKutta.hxx:60
void info(std::ostream &os)
Definition RungeKutta.hxx:490
The namespace for the Sandals library.
Definition Sandals.hh:73
Eigen::Vector< Real, Eigen::Dynamic > VectorX
Definition Sandals.hh:108
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Definition Sandals.hh:109
double Real
Definition Sandals.hh:84
static Real const EPSILON_HIGH
Definition Sandals.hh:123
static Real const SQRT_EPSILON
Definition Sandals.hh:121
int Integer
Definition Sandals.hh:85
Class container for the numerical solution of a system of ODEs/DAEs.
Definition Solution.hxx:57
Vector t
Definition Solution.hxx:62
void resize(Integer size)
Definition Solution.hxx:82
void conservative_resize(Integer size)
Definition Solution.hxx:92
MatrixN x
Definition Solution.hxx:63
Integer size() const
Definition Solution.hxx:121
MatrixM h
Definition Solution.hxx:64
Struct container for the Butcher tableau of a Runge-Kutta method.
Definition Tableau.hxx:36
bool is_embedded
Definition Tableau.hxx:49
enum class type :Integer {ERK=0, IRK=1, DIRK=2} Type
Definition Tableau.hxx:37
Eigen::Matrix< Real, S, S > Matrix
Definition Tableau.hxx:39
Integer order
Definition Tableau.hxx:43
Vector b_e
Definition Tableau.hxx:47
Type type
Definition Tableau.hxx:42
Eigen::Vector< Real, S > Vector
Definition Tableau.hxx:38
Matrix A
Definition Tableau.hxx:45
std::string name
Definition Tableau.hxx:41
Vector b
Definition Tableau.hxx:46
Vector c
Definition Tableau.hxx:48