Sandals  v0.0.0
A C++ library for ODEs/DAEs integration
Loading...
Searching...
No Matches
RungeKutta.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *\
2 * Copyright (c) 2025, Davide Stocco and Enrico Bertolazzi. *
3 * *
4 * The Sandals project is distributed under the BSD 2-Clause License. *
5 * *
6 * Davide Stocco Enrico Bertolazzi *
7 * University of Trento University of Trento *
8 * e-mail: davide.stocco@unitn.it e-mail: enrico.bertolazzi@unitn.it *
9\* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
10
11#pragma once
12
13#ifndef SANDALS_RUNGEKUTTA_HH
14#define SANDALS_RUNGEKUTTA_HH
15
16#include <Sandals.hh>
17
22
23#include <Sandals/Tableau.hh>
24#include <Sandals/Solution.hh>
25
26namespace Sandals {
27
28 /*\
29 | ____ _ __ _ _
30 | | _ \ _ _ _ __ __ _ ___| |/ / _| |_| |_ __ _
31 | | |_) | | | | '_ \ / _` |/ _ \ ' / | | | __| __/ _` |
32 | | _ <| |_| | | | | (_| | __/ . \ |_| | |_| || (_| |
33 | |_| \_\\__,_|_| |_|\__, |\___|_|\_\__,_|\__|\__\__,_|
34 | |___/
35 \*/
36
48 template <typename Real, Integer S, Integer N, Integer M>
50 {
51 using VectorX = Eigen::Vector<Real, Eigen::Dynamic>;
52 using MatrixX = Eigen::Matrix<Real, Eigen::Dynamic, Eigen::Dynamic>;
53 using VectorK = Eigen::Vector<Real, N*S>;
54 using MatrixK = Eigen::Matrix<Real, N, S>;
55 using MatrixJ = Eigen::Matrix<Real, N*S, N*S>;
56 using VectorP = Eigen::Matrix<Real, N+M, 1>;
57 using MatrixP = Eigen::Matrix<Real, N+M, N+M>;
58 using NewtonX = Optimist::RootFinder::Newton<Real, N>;
59 using NewtonK = Optimist::RootFinder::Newton<Real, N*S>;
66
67 public:
69 const Real SQRT_EPSILON{std::sqrt(EPSILON)}; \
70
72 using Type = typename Tableau<Real, S>::Type;
73 using Time = Eigen::Vector<Real, Eigen::Dynamic>;
74
75 private:
80 Real m_absolute_tolerance{EPSILON_HIGH};
81 Real m_relative_tolerance{EPSILON_HIGH};
82 Real m_safety_factor{0.9};
85 Real m_min_step{EPSILON_HIGH};
87 bool m_adaptive{true};
88 bool m_verbose{false};
89 bool m_reverse{false};
90
91 Eigen::FullPivLU<MatrixP> m_lu;
92 Real m_projection_tolerance{EPSILON_HIGH};
94 bool m_projection{true};
95
96 public:
100 RungeKutta(const RungeKutta &) = delete;
101
105 RungeKutta & operator=(RungeKutta const &) = delete;
106
111 RungeKutta(Tableau<Real, S> const &t_tableau)
112 : m_tableau(t_tableau) {
113 this->verbose_mode(this->m_verbose);
114 }
115
121 RungeKutta(Tableau<Real, S> const &t_tableau, System t_system)
122 : m_tableau(t_tableau), m_system(t_system) {
123 this->verbose_mode(this->m_verbose);
124 }
125
130 Type type() const {return this->m_tableau.type;}
131
136 bool is_erk() const {return this->m_tableau.type == Type::ERK;}
137
142 bool is_irk() const {return this->m_tableau.type == Type::IRK;}
143
148 bool is_dirk() const {return this->m_tableau.type == Type::DIRK;}
149
154 Tableau<Real, S> & tableau() {return this->m_tableau;}
155
160 Tableau<Real, S> const & tableau() const {return this->m_tableau;}
161
166 Integer stages() const {return S;}
167
172 std::string name() const {return this->m_tableau.name;}
173
178 Integer order() const {return this->m_tableau.order;}
179
184 bool is_embedded() const {return this->m_tableau.is_embedded;}
185
190 MatrixS A() const {return this->m_tableau.A;}
191
196 VectorS b() const {return this->m_tableau.b;}
197
202 VectorS b_embedded() const {return this->m_tableau.b_e;}
203
208 VectorS c() const {return c;}
209
214 System system() {return this->m_system;}
215
220 void system(System t_system) {this->m_system = t_system;}
221
238 ) {
239 this->m_system = std::make_shared<ImplicitWrapper<Real, N, M>>(
240 F, JF_x, JF_x_dot, h, Jh_x, in_domain);
241 }
242
254 std::string name,
261 ) {
262 this->m_system = std::make_shared<ImplicitWrapper<Real, N, M>>(
263 name, F, JF_x, JF_x_dot, h, Jh_x, in_domain);
264 }
265
280 ) {
281 this->m_system = std::make_shared<ExplicitWrapper<Real, N, M>>(
282 f, Jf_x, h, Jh_x, in_domain);
283 }
284
295 std::string name,
301 ) {
302 this->m_system = std::make_shared<ExplicitWrapper<Real, N, M>>(
303 name, f, Jf_x, h, Jh_x, in_domain);
304 }
305
322 ) {
323 this->m_system = std::make_shared<LinearWrapper<Real, N, M>>(
324 E, A, b, h, Jh_x, in_domain);
325 }
326
338 std::string name,
345 ) {
346 this->m_system = std::make_shared<LinearWrapper<Real, N, M>>(
347 name, E, A, b, h, Jh_x, in_domain);
348 }
349
372
397
402 bool has_system() {return this->m_system != nullptr;}
403
408 Real absolute_tolerance() {return this->m_absolute_tolerance;}
409
414 void absolute_tolerance(Real t_absolute_tolerance) {this->m_absolute_tolerance = t_absolute_tolerance;}
415
420 Real relative_tolerance() {return this->m_relative_tolerance;}
421
426 void relative_tolerance(Real t_relative_tolerance) {this->m_relative_tolerance = t_relative_tolerance;}
427
432 Real & safety_factor() {return this->m_safety_factor;}
433
438 void safety_factor(Real t_safety_factor) {this->m_safety_factor = t_safety_factor;}
439
444 Real & min_safety_factor() {return this->m_min_safety_factor;}
445
450 void min_safety_factor(Real t_min_safety_factor) {this->m_min_safety_factor = t_min_safety_factor;}
451
456 Real & max_safety_factor() {return this->m_max_safety_factor;}
457
462 void max_safety_factor(Real t_max_safety_factor) {this->m_max_safety_factor = t_max_safety_factor;}
463
468 Real & min_step() {return this->m_min_step;}
469
474 void min_step(Real t_min_step) {this->m_min_step = t_min_step;}
475
480 Integer & max_substeps() {return this->m_max_substeps;}
481
486 void max_substeps(Integer t_max_substeps) {this->m_max_substeps = t_max_substeps;}
487
492 bool adaptive_mode() {return this->m_adaptive;}
493
498 void adaptive(bool t_adaptive) {this->m_adaptive = t_adaptive;}
499
503 void enable_adaptive_mode() {this->m_adaptive = true;}
504
508 void disable_adaptive_mode() {this->m_adaptive = false;}
509
514 bool verbose_mode() {return this->m_verbose;}
515
520 void verbose_mode(bool t_verbose) {
521 this->m_verbose = t_verbose;
522 this->m_newtonX.verbose_mode(t_verbose);
523 this->m_newtonK.verbose_mode(t_verbose);
524 }
525
529 void enable_verbose_mode() {this->verbose_mode(true);}
530
534 void disable_verbose_mode() {this->verbose_mode(false);}
535
540 bool reverse_mode() {return this->m_reverse;}
541
546 void reverse(bool t_reverse) {this->m_reverse = t_reverse;}
547
551 void enable_reverse_mode() {this->m_reverse = true;}
552
556 void disable_reverse_mode() {this->m_reverse = false;}
557
562 Real projection_tolerance() {return this->m_projection_tolerance;}
563
568 void projection_tolerance(Real t_projection_tolerance)
569 {this->m_projection_tolerance = t_projection_tolerance;}
570
575 Integer & max_projection_iterations() {return this->m_max_projection_iterations;}
576
581 void max_projection_iterations(Integer t_max_projection_iterations)
582 {this->m_max_projection_iterations = t_max_projection_iterations;}
583
588 bool projection() {return this->m_projection;}
589
594 void projection(bool t_projection) {this->m_projection = t_projection;}
595
599 void enable_projection() {this->m_projection = true;}
600
604 void disable_projection() {this->m_projection = false;}
605
640 Real estimate_step(VectorN const &x, VectorN const &x_e, Real h_k) const
641 {
642 Real desired_error{this->m_absolute_tolerance + this->m_relative_tolerance *
643 std::max(x.array().abs().maxCoeff(), x_e.array().abs().maxCoeff())};
644 Real truncation_error{(x - x_e).array().abs().maxCoeff()};
645 return h_k * std::min(this->m_max_safety_factor, std::max(this->m_min_safety_factor,
646 this->m_safety_factor * std::pow(desired_error/truncation_error,
647 1.0/std::max(this->m_tableau.order, this->m_tableau.order_e))));
648 }
649
654 std::string info() const {
655 std::ostringstream os;
656 os
657 << "Runge-Kutta method:\t" << this->name() << std::endl
658 << "\t- order:\t" << this->order() << std::endl
659 << "\t- stages:\t" << this->stages() << std::endl
660 << "\t- type:\t";
661 switch (this->type()) {
662 case Type::ERK: os << "explicit"; break;
663 case Type::IRK: os << "implicit"; break;
664 case Type::DIRK: os << "diagonally implicit"; break;
665 }
666 os
667 << std::endl
668 << "\t- embedded:\t" << this->is_embedded() << std::endl;
669 if (this->has_system()) {
670 os << "\t- system:\t" << this->m_system->name() << std::endl;
671 } else {
672 os << "\t- system:\t" << "none" << std::endl;
673 }
674 return os.str();
675 }
676
681 void info(std::ostream &os) {os << this->info();}
682
683 /*\
684 | _____ ____ _ __
685 | | ____| _ \| |/ /
686 | | _| | |_) | ' /
687 | | |___| _ <| . \
688 | |_____|_| \_\_|\_\
689 |
690 \*/
691
705 bool erk_explicit_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new) const
706 {
707 using Eigen::all;
708 using Eigen::seqN;
709
710 // Compute the K variables in the case of an explicit method and explicit system
711 VectorN x_node;
712 MatrixK K;
713 for (Integer i{0}; i < S; ++i) {
714 x_node = x_old + K(all, seqN(0, i)) * this->m_tableau.A(i, seqN(0, i)).transpose();
715 if (!this->m_reverse) {
716 K.col(i) = h_old * static_cast<Explicit<Real, N, M> const *>(this->m_system.get())->f(x_node, t_old + h_old*this->m_tableau.c(i));
717 } else {
718 K.col(i) = h_old * static_cast<Explicit<Real, N, M> const *>(this->m_system.get())->f_reverse(x_node, t_old + h_old*this->m_tableau.c(i));
719 }
720 }
721 if (!K.allFinite()) {return false;}
722
723 // Perform the step and obtain the next state
724 x_new = x_old + K * this->m_tableau.b;
725
726 // Adapt next step
727 if (this->m_adaptive && this->m_tableau.is_embedded) {
728 VectorN x_emb = x_old + K * this->m_tableau.b_e;
729 h_new = this->estimate_step(x_new, x_emb, h_old);
730 }
731 return true;
732 }
733
754 void erk_implicit_function(Integer s, VectorN const &x, Real t, Real h, MatrixK const &K, VectorN &fun) const
755 {
756 using Eigen::all;
757 using Eigen::seqN;
758 VectorN x_node(x + K(all, seqN(0, s)) * this->m_tableau.A(s, seqN(0, s)).transpose());
759 if (!this->m_reverse) {
760 fun = this->m_system->F(x_node, K.col(s)/h, t + h * this->m_tableau.c(s));
761 } else {
762 fun = this->m_system->F_reverse(x_node, K.col(s)/h, t + h * this->m_tableau.c(s));
763 }
764 }
765
797 void erk_implicit_jacobian(Integer s, VectorN const &x, Real t, Real h, MatrixK const &K, MatrixN &jac) const
798 {
799 using Eigen::all;
800 using Eigen::seqN;
801 VectorN x_node(x + K(all, seqN(0, s)) * this->m_tableau.A(s, seqN(0, s)).transpose());
802 if (!this->m_reverse) {
803 jac = this->m_system->JF_x_dot(x_node, K.col(s)/h, t + h * this->m_tableau.c(s)) / h;
804 } else {
805 jac = this->m_system->JF_x_dot_reverse(x_node, K.col(s)/h, t + h * this->m_tableau.c(s)) / h;
806 }
807 }
808
822 bool erk_implicit_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
823 {
824 MatrixK K;
825 VectorN K_sol;
826 VectorN K_ini(VectorN::Zero());
827
828 // Check if the solver converged
829 for (Integer s{0}; s < S; ++s) {
830 if (this->m_newtonX.solve(
831 [this, s, &K, &x_old, t_old, h_old](VectorN const &K_fun, VectorN &fun)
832 {K.col(s) = K_fun; this->erk_implicit_function(s, x_old, t_old, h_old, K, fun);},
833 [this, s, &K, &x_old, t_old, h_old](VectorN const &K_jac, MatrixN &jac)
834 {K.col(s) = K_jac; this->erk_implicit_jacobian(s, x_old, t_old, h_old, K, jac);},
835 K_ini, K_sol)) {
836 K.col(s) = K_sol;
837 } else {
838 return false;
839 }
840 }
841
842 // Perform the step and obtain the next state
843 x_new = x_old + K * this->m_tableau.b;
844
845 // Adapt next step
846 if (this->m_adaptive && this->m_tableau.is_embedded) {
847 VectorN x_emb(x_old + K * this->m_tableau.b_e);
848 h_new = this->estimate_step(x_new, x_emb, h_old);
849 }
850 return true;
851 }
852
853 /*\
854 | ___ ____ _ __
855 | |_ _| _ \| |/ /
856 | | || |_) | ' /
857 | | || _ <| . \
858 | |___|_| \_\_|\_\
859 |
860 \*/
861
886 void irk_function(VectorN const &x, Real t, Real h, VectorK const &K, VectorK &fun) const
887 {
888 VectorN x_node;
889 MatrixK K_mat{K.reshaped(N, S)};
890 MatrixK fun_mat;
891 for (Integer i{0}; i < S; ++i) {
892 x_node = x + K_mat * this->m_tableau.A.row(i).transpose();
893 if (!this->m_reverse) {
894 fun_mat.col(i) = this->m_system->F(x_node, K_mat.col(i)/h, t + h * this->m_tableau.c(i));
895 } else {
896 fun_mat.col(i) = this->m_system->F_reverse(x_node, K_mat.col(i)/h, t + h * this->m_tableau.c(i));
897 }
898 }
899 fun = fun_mat.reshaped(N*S, 1);
900 }
901
938 void irk_jacobian(VectorN const &x, Real t, Real h, VectorK const &K, MatrixJ &jac) const
939 {
940 using Eigen::seqN;
941
942 // Reset the Jacobian matrix
943 jac.setZero();
944
945 // Loop through each equation of the system
946 MatrixK K_mat{K.reshaped(N, S)};
947 Real t_node;
948 VectorN x_node, x_dot_node;
949 MatrixN JF_x, JF_x_dot;
950 auto idx = seqN(0, N), jdx = seqN(0, N);
951 for (Integer i{0}; i < S; ++i) {
952 t_node = t + h * this->m_tableau.c(i);
953 x_node = x + K_mat * this->m_tableau.A.row(i).transpose();
954
955 // Compute the Jacobians with respect to x and x_dot
956 x_dot_node = K_mat.col(i) / h;
957 if (!this->m_reverse) {
958 JF_x = this->m_system->JF_x(x_node, x_dot_node, t_node);
959 JF_x_dot = this->m_system->JF_x_dot(x_node, x_dot_node, t_node);
960 } else {
961 JF_x = this->m_system->JF_x_reverse(x_node, x_dot_node, t_node);
962 JF_x_dot = this->m_system->JF_x_dot_reverse(x_node, x_dot_node, t_node);
963 }
964
965 // Combine the Jacobians with respect to x and x_dot to obtain the Jacobian with respect to K
966 idx = seqN(i*N, N);
967 for (Integer j{0}; j < S; ++j) {
968 jdx = seqN(j*N, N);
969 if (i == j) {
970 jac(idx, jdx) = this->m_tableau.A(i,j) * JF_x + JF_x_dot / h;
971 } else {
972 jac(idx, jdx) = this->m_tableau.A(i,j) * JF_x;
973 }
974 }
975 }
976 }
977
991 bool irk_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
992 {
993 VectorK K;
994 VectorK K_ini(VectorK::Zero());
995
996 // Check if the solver converged
997 if (!this->m_newtonK.solve(
998 [this, &x_old, t_old, h_old](VectorK const &K_fun, VectorK &fun)
999 {this->irk_function(x_old, t_old, h_old, K_fun, fun);},
1000 [this, &x_old, t_old, h_old](VectorK const &K_jac, MatrixJ &jac)
1001 {this->irk_jacobian(x_old, t_old, h_old, K_jac, jac);},
1002 K_ini, K))
1003 {return false;}
1004
1005 // Perform the step and obtain the next state
1006 x_new = x_old + K.reshaped(N, S) * this->m_tableau.b;
1007
1008 // Adapt next step
1009 if (this->m_adaptive && this->m_tableau.is_embedded) {
1010 VectorN x_emb(x_old + K.reshaped(N, S) * this->m_tableau.b_e);
1011 h_new = this->estimate_step(x_new, x_emb, h_old);
1012 }
1013 return true;
1014 }
1015
1016 /*\
1017 | ____ ___ ____ _ __
1018 | | _ \_ _| _ \| |/ /
1019 | | | | | || |_) | ' /
1020 | | |_| | || _ <| . \
1021 | |____/___|_| \_\_|\_\
1022 |
1023 \*/
1024
1045 void dirk_function(Integer n, VectorN const &x, Real t, Real h, MatrixK const &K, VectorN &fun) const
1046 {
1047 using Eigen::all;
1048 using Eigen::seqN;
1049 VectorN x_node(x + K(all, seqN(0, n+1)) * this->m_tableau.A(n, seqN(0, n+1)).transpose());
1050 if (!this->m_reverse) {
1051 fun = this->m_system->F(x_node, K.col(n)/h, t + h * this->m_tableau.c(n));
1052 } else {
1053 fun = this->m_system->F_reverse(x_node, K.col(n)/h, t + h * this->m_tableau.c(n));
1054 }
1055 }
1056
1090 void dirk_jacobian(Integer n, VectorN const &x, Real t, Real h, MatrixK const &K, MatrixN &jac) const
1091 {
1092 using Eigen::all;
1093 using Eigen::seqN;
1094 Real t_node{t + h * this->m_tableau.c(n)};
1095 VectorN x_node(x + K(all, seqN(0, n+1)) * this->m_tableau.A(n, seqN(0, n+1)).transpose());
1096 VectorN x_dot_node(K.col(n)/h);
1097 if (!this->m_reverse) {
1098 jac = this->m_tableau.A(n,n) * this->m_system->JF_x(x_node, x_dot_node, t_node) +
1099 this->m_system->JF_x_dot(x_node, x_dot_node, t_node) / h;
1100 } else {
1101 jac = this->m_tableau.A(n,n) * this->m_system->JF_x_reverse(x_node, x_dot_node, t_node) +
1102 this->m_system->JF_x_dot_reverse(x_node, x_dot_node, t_node) / h;
1103 }
1104 }
1105
1119 bool dirk_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
1120 {
1121 MatrixK K;
1122 VectorN K_sol;
1123 VectorN K_ini(VectorN::Zero());
1124
1125 // Check if the solver converged at each step
1126 for (Integer n{0}; n < S; ++n) {
1127 if (this->m_newtonX.solve(
1128 [this, n, &K, &x_old, t_old, h_old](VectorN const &K_fun, VectorN &fun)
1129 {K.col(n) = K_fun; this->dirk_function(n, x_old, t_old, h_old, K, fun);},
1130 [this, n, &K, &x_old, t_old, h_old](VectorN const &K_jac, MatrixN &jac)
1131 {K.col(n) = K_jac; this->dirk_jacobian(n, x_old, t_old, h_old, K, jac);},
1132 K_ini, K_sol)) {
1133 K.col(n) = K_sol;
1134 } else {
1135 return false;
1136 }
1137 }
1138
1139 // Perform the step and obtain the next state
1140 x_new = x_old + K * this->m_tableau.b;
1141
1142 // Adapt next step
1143 if (this->m_adaptive && this->m_tableau.is_embedded) {
1144 VectorN x_emb(x_old + K * this->m_tableau.b_e);
1145 h_new = this->estimate_step(x_new, x_emb, h_old);
1146 }
1147 return true;
1148 }
1149
1161 bool step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
1162 {
1163 #define CMD "Sandals::RungeKutta::step(...): "
1164
1165 SANDALS_ASSERT(this->m_system->in_domain(x_old, t_old), CMD "in " << this->m_tableau.name <<
1166 " solver, at t = " << t_old << ", x = " << x_old.transpose() << ", system out of domain.");
1167
1168 if (this->is_erk() && this->m_system->is_explicit()) {
1169 return this->erk_explicit_step(x_old, t_old, h_old, x_new, h_new);
1170 } else if (this->is_erk() && this->m_system->is_implicit()) {
1171 return this->erk_implicit_step(x_old, t_old, h_old, x_new, h_new);
1172 } else if (this->is_dirk()) {
1173 return this->dirk_step(x_old, t_old, h_old, x_new, h_new);
1174 } else {
1175 return this->irk_step(x_old, t_old, h_old, x_new, h_new);
1176 }
1177
1178 #undef CMD
1179 }
1180
1194 bool advance(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
1195 {
1196 #define CMD "Sandals::RungeKutta::advance(...): "
1197
1198 // Check step size
1199 SANDALS_ASSERT(h_old > Real(0.0), CMD "in " << this->m_tableau.name << " solver, h = "<<
1200 h_old << ", expected > 0.");
1201
1202 // If the integration step failed, try again with substepping
1203 if (!this->step(x_old, t_old, h_old, x_new, h_new))
1204 {
1205 VectorN x_tmp(x_old);
1206 Real t_tmp{t_old}, h_tmp{h_old / Real(2.0)};
1207
1208 // Substepping logic
1209 Integer max_k{this->m_max_substeps * this->m_max_substeps}, k{2};
1210 Real h_new_tmp;
1211 while (k > 0) {
1212 // Calculate the next step with substepping logic
1213 if (this->step(x_tmp, t_tmp, h_tmp, x_new, h_new_tmp)) {
1214
1215 // Accept the step
1216 h_tmp = h_new_tmp;
1217
1218 // If substepping is enabled, double the step size
1219 if (k > 0 && k < max_k) {
1220 k -= 1;
1221 // If the substepping index is even, double the step size
1222 if (k % 2 == 0) {
1223 h_tmp = Real(2.0) * h_tmp;
1224 if (this->m_verbose) {
1225 SANDALS_WARNING(CMD "in " << this->m_tableau.name << " solver, at t = " << t_tmp <<
1226 ", integration succedded disable one substepping layer.");
1227 }
1228 }
1229 }
1230
1231 // Check the infinity norm of the solution
1232 SANDALS_ASSERT(std::isfinite(x_tmp.maxCoeff()), CMD "in " << this->m_tableau.name <<
1233 " solver, at t = " << t_tmp << ", ||x||_inf = inf, computation interrupted.");
1234
1235 } else {
1236
1237 // If the substepping index is too high, abort the integration
1238 k += 2;
1239 SANDALS_ASSERT(k < max_k, CMD "in " << this->m_tableau.name << " solver, at t = " <<
1240 t_tmp << ", integration failed with h = " << h_tmp << ", aborting.");
1241 return false;
1242
1243 // Otherwise, try again with a smaller step
1244 if (this->m_verbose) {SANDALS_WARNING(CMD "in " << this->m_tableau.name << " solver, " <<
1245 "at t = " << t_tmp << ", integration failed, adding substepping layer.");}
1246 h_tmp /= Real(2.0);
1247 continue;
1248 }
1249
1250 // Store independent variable (or time)
1251 t_tmp += h_tmp;
1252 }
1253
1254 // Store output states substepping solutions
1255 x_new = x_tmp;
1256 h_new = h_tmp;
1257 }
1258
1259 // Project intermediate solution on the invariants
1260 if (this->m_projection) {
1261 VectorN x_projected;
1262 if (this->project(x_new, t_old + h_new, x_projected)) {
1263 x_new = x_projected;
1264 return true;
1265 } else {
1266 return false;
1267 }
1268 } else {
1269 return true;
1270 }
1271
1272 #undef CMD
1273 }
1274
1284 bool solve(VectorX const &t_mesh, VectorN const &ics, Solution<Real, N, M> &sol)
1285 {
1286 using Eigen::last;
1287
1288 // Instantiate output
1289 sol.resize(t_mesh.size());
1290
1291 // Store first step
1292 sol.t(0) = t_mesh(0);
1293 sol.x.col(0) = ics;
1294 sol.h.col(0) = this->m_system->h(ics, t_mesh(0));
1295
1296 // Update the current step
1297 Integer step{0};
1298 VectorN x_step{ics};
1299 Real t_step{t_mesh(0)}, h_step{t_mesh(1) - t_mesh(0)}, h_tmp_step{h_step}, h_new_step;
1300 bool mesh_point_bool, saturation_bool;
1301
1302 while (true) {
1303 // Integrate system
1304 if (!this->advance(sol.x.col(step), t_step, h_step, x_step, h_new_step)) {return false;}
1305
1306 // Update the current step
1307 t_step += h_step;
1308
1309 // Saturate the suggested timestep
1310 mesh_point_bool = std::abs(t_step - t_mesh(step+1)) < SQRT_EPSILON;
1311 saturation_bool = t_step + h_new_step > t_mesh(step+1) + SQRT_EPSILON;
1312 if (this->m_adaptive && this->m_tableau.is_embedded && !mesh_point_bool && saturation_bool) {
1313 h_tmp_step = h_new_step; // Used to store the previous step and keep the integration pace
1314 h_step = t_mesh(step+1) - t_step;
1315 } else {
1316 h_step = h_new_step;
1317 }
1318
1319 // Store solution if the step is a mesh point
1320 if (!this->m_adaptive || mesh_point_bool) {
1321 // Update temporaries
1322 step += 1;
1323 h_step = h_tmp_step;
1324
1325 // Update outputs
1326 sol.t(step) = t_step;
1327 sol.x.col(step) = x_step;
1328 sol.h.col(step) = this->m_system->h(x_step, t_step);
1329
1330 // Check if the current step is the last one
1331 if (std::abs(t_step - t_mesh(last)) < SQRT_EPSILON) {break;}
1332 }
1333 }
1334 return true;
1335 }
1336
1346 bool adaptive_solve(VectorX const &t_mesh, VectorN const &ics, Solution<Real, N, M> &sol)
1347 {
1348 using Eigen::all;
1349 using Eigen::last;
1350
1351 #define CMD "Sandals::RungeKutta::adaptive_solve(...): "
1352
1353 // Check if the adaptive method is enabled and the method is embedded
1354 if (!this->is_embedded()) {
1355 SANDALS_WARNING(CMD "the method is not embedded, using solve(...) method.");
1356 return this->solve(t_mesh, ics, sol);
1357 } else if (!this->m_adaptive) {
1358 SANDALS_WARNING(CMD "adaptive method is disabled, using solve(...) method.");
1359 return this->solve(t_mesh, ics, sol);
1360 }
1361
1362 // Instantiate output
1363 Real h_step{t_mesh(1) - t_mesh(0)}, h_new_step, scale{100.0};
1364 Real h_min{std::max(this->m_min_step, h_step/scale)}, h_max{scale*h_step};
1365 if (this->m_tableau.is_embedded) {
1366 Integer safety_length{static_cast<Integer>(std::ceil(std::abs(t_mesh(last) - t_mesh(0))/(2.0*h_min)))};
1367 sol.resize(safety_length);
1368 } else {
1369 sol.resize(t_mesh.size());
1370 }
1371
1372 // Store first step
1373 sol.t(0) = t_mesh(0);
1374 sol.x.col(0) = ics;
1375 sol.h.col(0) = this->m_system->h(ics, t_mesh(0));
1376
1377 // Instantiate temporary variables
1378 Integer step{0};
1379 VectorN x_step{ics};
1380
1381 while (true) {
1382 // Integrate system
1383 this->advance(sol.x.col(step), sol.t(step), h_step, x_step, h_new_step);
1384
1385 // Saturate the suggested timestep
1386 if (this->m_adaptive && this->m_tableau.is_embedded) {
1387 h_step = std::max(std::min(h_new_step, h_max), h_min);
1388 }
1389
1390 SANDALS_ASSERT(step < sol.size(), CMD "safety length exceeded.");
1391
1392 // Store solution
1393 sol.t(step+1) = sol.t(step) + h_step;
1394 sol.x.col(step+1) = x_step;
1395 sol.h.col(step+1) = this->m_system->h(x_step, sol.t(step+1));
1396
1397 // Check if the current step is the last one
1398 if (sol.t(step+1) + h_step > t_mesh(last)) {break;}
1399
1400 // Update steps counter
1401 step += 1;
1402 }
1403
1404 // Resize the output
1406
1407 return true;
1408
1409 #undef CMD
1410 }
1411
1421 bool project(VectorN const &x, Real t, VectorN &x_projected)
1422 {
1423 #define CMD "Sandals::RungeKutta::project(...): "
1424
1425 // Check if there are any constraints
1426 x_projected = x;
1427 if (M > Integer(0)) {
1428
1429 VectorM h;
1430 MatrixM Jh_x;
1431 VectorP b, x_step;
1432 MatrixP A;
1433 A.setZero();
1434 A.template block<N, N>(0, 0) = MatrixN::Identity();
1435 for (Integer k{0}; k < this->m_max_projection_iterations; ++k) {
1436
1437 /* Standard projection method
1438 [A] {x} = {b}
1439 / I Jh_x^T \ / dx \ / x_t - x_k \
1440 | | | | = | |
1441 \ Jh_x 0 / \ lambda / \ -h /
1442 */
1443
1444 // Evaluate the invariants vector and its Jacobian
1445 h = this->m_system->h(x_projected, t);
1446 Jh_x = this->m_system->Jh_x(x_projected, t);
1447
1448 // Check if the solution is found
1449 if (h.norm() < this->m_projection_tolerance) {return true;}
1450
1451 // Build the linear system
1452 A.template block<N, M>(0, N) = Jh_x.transpose();
1453 A.template block<M, N>(N, 0) = Jh_x;
1454 b.template head<N>() = x - x_projected;
1455 b.template tail<M>() = -h;
1456
1457 // Compute the solution of the linear system
1458 this->m_lu.compute(A);
1459 SANDALS_ASSERT(this->m_lu.rank() == N+M, CMD "singular Jacobian detected.");
1460 x_step = this->m_lu.solve(b);
1461
1462 // Check if the step is too small
1463 if (x_step.norm() < this->m_projection_tolerance * this->m_projection_tolerance) {return false;}
1464
1465 // Update the solution
1466 x_projected.noalias() += x_step(Eigen::seqN(0, N));
1467 }
1468 if (this->m_verbose) {SANDALS_WARNING(CMD "maximum number of iterations reached.");}
1469 return false;
1470 } else {
1471 return true;
1472 }
1473
1474 #undef CMD
1475 }
1476
1488 bool project_ics(VectorN const &x, Real t, std::vector<Integer> const & projected_equations,
1489 std::vector<Integer> const & projected_invariants, VectorN &x_projected) const
1490 {
1491 #define CMD "Sandals::RungeKutta::project_ics(...): "
1492
1493 Integer X{static_cast<Integer>(projected_equations.size())};
1494 Integer H{static_cast<Integer>(projected_invariants.size())};
1495
1496 // Check if there are any constraints
1497 x_projected = x;
1498 if (H > Integer(0)) {
1499
1500 VectorM h;
1501 MatrixM Jh_x;
1502 VectorX b(X+H), x_step(X+H);
1503 MatrixX A(X+H, X+H);
1504 A.setZero();
1505 A.block(0, 0, X, X) = MatrixX::Identity(X+H, X+H);
1506 Eigen::FullPivLU<MatrixX> lu;
1507 for (Integer k{0}; k < this->m_max_projection_iterations; ++k) {
1508
1509 /* Standard projection method
1510 [A] {x} = {b}
1511 / I Jh_x^T \ / dx \ / x_t - x_k \
1512 | | | | = | |
1513 \ Jh_x 0 / \ lambda / \ -h /
1514 */
1515
1516 // Evaluate the invariants vector and its Jacobian
1517 h = this->m_system->h(x_projected, t);
1518 Jh_x = this->m_system->Jh_x(x_projected, t);
1519
1520 // Select only the projected invariants
1521 h = h(projected_invariants);
1522 Jh_x = Jh_x(projected_invariants, projected_equations);
1523
1524 // Check if the solution is found
1525 if (h.norm() < this->m_projection_tolerance) {return true;}
1526
1527 // Build the linear system
1528 A.block(0, X, X, H) = Jh_x.transpose();
1529 A.block(X, 0, H, X) = Jh_x;
1530 b.head(X) = x(projected_equations) - x_projected(projected_equations);
1531 b.tail(H) = -h;
1532
1533 // Compute the solution of the linear system
1534 lu.compute(A);
1535 SANDALS_ASSERT(lu.rank() == X+H, CMD "singular Jacobian detected.");
1536 x_step = this->m_lu.solve(b);
1537
1538 // Check if the step is too small
1539 if (x_step.norm() < this->m_projection_tolerance * this->m_projection_tolerance) {return false;}
1540
1541 // Update the solution
1542 x_projected(projected_equations).noalias() += x_step;
1543 }
1544 if (this->m_verbose) {SANDALS_WARNING(CMD "maximum number of iterations reached.");}
1545 return false;
1546 } else {
1547 return true;
1548 }
1549
1550 #undef CMD
1551 }
1552
1560 Real estimate_order(std::vector<VectorX> const &t_mesh, VectorN const &ics, std::function<MatrixX(VectorX)> &sol)
1561 {
1562 using Eigen::last;
1563
1564 #define CMD "Sandals::RungeKutta::estimate_order(...): "
1565
1566 SANDALS_ASSERT(t_mesh.size() > Integer(1), CMD "expected at least two time meshes.");
1567
1568 for (Integer i{0}; i < static_cast<Integer>(t_mesh.size()); ++i) {
1569
1570 // Check if the time mesh has the same initial and final time
1571 SANDALS_ASSERT((t_mesh[0](0) - t_mesh[i](0)) < SQRT_EPSILON,
1572 CMD "expected the same initial time.");
1573 SANDALS_ASSERT((t_mesh[0](last) - t_mesh[i](last)) < SQRT_EPSILON,
1574 CMD "expected the same final time.");
1575
1576 // Check if the mesh step is fixed all over the time mesh
1577 for (Integer j{1}; j < static_cast<Integer>(t_mesh[i].size()); ++j) {
1578 SANDALS_ASSERT((t_mesh[i](j) - t_mesh[i](j-1)) - (t_mesh[i](1) - t_mesh[i](0)) < SQRT_EPSILON,
1579 CMD "expected a fixed step.");
1580 }
1581 }
1582
1583 // Solve the system for each time scale
1584 Solution<Real, N, M> sol_num;
1585 MatrixX sol_ana;
1586 VectorX h_vec(t_mesh.size()), e_vec(t_mesh.size());
1587 for (Integer i{0}; i < static_cast<Integer>(t_mesh.size()); ++i) {
1588 SANDALS_ASSERT(this->solve(t_mesh[i], ics, sol_num), CMD "failed to solve the system for " <<
1589 "the" << i << "-th time mesh.");
1590 sol_ana = sol(sol_num.t);
1591 SANDALS_ASSERT(sol_ana.rows() == sol_num.x.rows(),
1592 CMD "expected the same number of states in analytical solution.");
1593 SANDALS_ASSERT(sol_ana.cols() == sol_num.x.cols(),
1594 CMD "expected the same number of steps in analytical solution.");
1595 h_vec(i) = std::abs(sol_num.t(1) - sol_num.t(0));
1596 e_vec(i) = (sol_ana - sol_num.x).array().abs().maxCoeff();
1597 }
1598
1599 // Compute the order of the method thorugh least squares
1600 VectorX A(h_vec.array().log());
1601 VectorX b(e_vec.array().log());
1602 return ((A.transpose() * A).ldlt().solve(A.transpose() * b))(0);
1603
1604 #undef CMD
1605 }
1606
1607 }; // class RungeKutta
1608
1609} // namespace Sandals
1610
1611#endif // SANDALS_RUNGEKUTTA_HH
#define CMD
#define SANDALS_BASIC_CONSTANTS(Real)
Definition Sandals.hh:70
#define SANDALS_ASSERT(COND, MSG)
Definition Sandals.hh:44
#define SANDALS_WARNING(MSG)
Definition Sandals.hh:53
Class container for the system of explicit ODEs.
Definition Explicit.hh:42
VectorF f_reverse(VectorF const &x, Real t) const
Definition Explicit.hh:160
virtual VectorF f(VectorF const &x, Real t) const =0
static const FunctionH DefaultH
Definition Explicit.hh:258
std::function< MatrixJH(VectorF const &, Real)> FunctionJH
Definition Explicit.hh:255
static const FunctionID DefaultID
Definition Explicit.hh:260
std::function< VectorF(VectorF const &, Real)> FunctionF
Definition Explicit.hh:252
std::function< MatrixJF(VectorF const &, Real)> FunctionJF
Definition Explicit.hh:253
std::function< VectorH(VectorF const &, Real)> FunctionH
Definition Explicit.hh:254
static const FunctionJH DefaultJH
Definition Explicit.hh:259
std::function< bool(VectorF const &, Real)> FunctionID
Definition Explicit.hh:256
Eigen::Matrix< Real, N, N > MatrixJF
Definition Implicit.hh:49
Eigen::Vector< Real, N > VectorF
Definition Implicit.hh:48
Eigen::Vector< Real, M > VectorH
Definition Implicit.hh:50
std::shared_ptr< Implicit< Real, N, M > > Pointer
Definition Implicit.hh:47
Eigen::Matrix< Real, M, N > MatrixJH
Definition Implicit.hh:51
std::function< bool(VectorF const &, Real)> FunctionID
Definition Implicit.hh:283
std::function< VectorH(VectorF const &, Real)> FunctionH
Definition Implicit.hh:281
std::function< VectorF(VectorF const &, VectorF const &, Real)> FunctionF
Definition Implicit.hh:279
std::function< MatrixJF(VectorF const &, VectorF const &, Real)> FunctionJF
Definition Implicit.hh:280
static const FunctionID DefaultID
Definition Implicit.hh:287
static const FunctionH DefaultH
Definition Implicit.hh:285
std::function< MatrixJH(VectorF const &, Real)> FunctionJH
Definition Implicit.hh:282
static const FunctionJH DefaultJH
Definition Implicit.hh:286
std::function< MatrixJH(VectorF const &, Real)> FunctionJH
Definition Linear.hh:219
std::function< VectorB(Real)> FunctionB
Definition Linear.hh:217
static const FunctionH DefaultH
Definition Linear.hh:222
std::function< bool(VectorF const &, Real)> FunctionID
Definition Linear.hh:220
static const FunctionJH DefaultJH
Definition Linear.hh:223
std::function< MatrixA(Real)> FunctionA
Definition Linear.hh:216
std::function< VectorH(VectorF const &, Real)> FunctionH
Definition Linear.hh:218
static const FunctionID DefaultID
Definition Linear.hh:224
std::function< MatrixE(Real)> FunctionE
Definition Linear.hh:215
Real & min_step()
Definition RungeKutta.hh:468
bool projection()
Definition RungeKutta.hh:588
bool adaptive_mode()
Definition RungeKutta.hh:492
bool is_erk() const
Definition RungeKutta.hh:136
typename Tableau< Real, S >::Matrix MatrixS
Definition RungeKutta.hh:61
void enable_reverse_mode()
Definition RungeKutta.hh:551
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.hh:1488
void projection(bool t_projection)
Definition RungeKutta.hh:594
typename Tableau< Real, S >::Vector VectorS
Definition RungeKutta.hh:60
Real & min_safety_factor()
Definition RungeKutta.hh:444
void semi_explicit_system(typename SemiExplicitWrapper< Real, N, M >::FunctionA A, typename SemiExplicitWrapper< Real, N, M >::FunctionTA TA_x, typename SemiExplicitWrapper< Real, N, M >::FunctionB b, typename SemiExplicitWrapper< Real, N, M >::FunctionJB Jb_x, typename SemiExplicitWrapper< Real, N, M >::FunctionH h=SemiExplicitWrapper< Real, N, M >::DefaultH, typename SemiExplicitWrapper< Real, N, M >::FunctionJH Jh_x=SemiExplicitWrapper< Real, N, M >::DefaultJH, typename SemiExplicitWrapper< Real, N, M >::FunctionID in_domain=SemiExplicitWrapper< Real, N, M >::DefaultID)
Definition RungeKutta.hh:360
bool reverse_mode()
Definition RungeKutta.hh:540
typename Implicit< Real, N, M >::MatrixJF MatrixN
Definition RungeKutta.hh:63
Real & safety_factor()
Definition RungeKutta.hh:432
Real m_absolute_tolerance
Definition RungeKutta.hh:80
Real m_min_safety_factor
Definition RungeKutta.hh:83
RungeKutta & operator=(RungeKutta const &)=delete
void irk_function(VectorN const &x, Real t, Real h, VectorK const &K, VectorK &fun) const
Definition RungeKutta.hh:886
typename Implicit< Real, N, M >::MatrixJH MatrixM
Definition RungeKutta.hh:65
Eigen::Matrix< Real, N+M, 1 > VectorP
Definition RungeKutta.hh:56
Tableau< Real, S > & tableau()
Definition RungeKutta.hh:154
System system()
Definition RungeKutta.hh:214
VectorS b_embedded() const
Definition RungeKutta.hh:202
RungeKutta(const RungeKutta &)=delete
bool step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hh:1161
typename Implicit< Real, N, M >::Pointer System
Definition RungeKutta.hh:71
void enable_projection()
Definition RungeKutta.hh:599
RungeKutta(Tableau< Real, S > const &t_tableau, System t_system)
Definition RungeKutta.hh:121
Real m_relative_tolerance
Definition RungeKutta.hh:81
void max_safety_factor(Real t_max_safety_factor)
Definition RungeKutta.hh:462
Integer & max_substeps()
Definition RungeKutta.hh:480
typename Implicit< Real, N, M >::VectorF VectorN
Definition RungeKutta.hh:62
void linear_system(std::string name, typename LinearWrapper< Real, N, M >::FunctionE E, typename LinearWrapper< Real, N, M >::FunctionA A, typename LinearWrapper< Real, N, M >::FunctionB b, typename LinearWrapper< Real, N, M >::FunctionH h=LinearWrapper< Real, N, M >::DefaultH, typename LinearWrapper< Real, N, M >::FunctionJH Jh_x=LinearWrapper< Real, N, M >::DefaultJH, typename LinearWrapper< Real, N, M >::FunctionID in_domain=LinearWrapper< Real, N, M >::DefaultID)
Definition RungeKutta.hh:337
RungeKutta(Tableau< Real, S > const &t_tableau)
Definition RungeKutta.hh:111
Eigen::Matrix< Real, N+M, N+M > MatrixP
Definition RungeKutta.hh:57
void adaptive(bool t_adaptive)
Definition RungeKutta.hh:498
Integer order() const
Definition RungeKutta.hh:178
void dirk_function(Integer n, VectorN const &x, Real t, Real h, MatrixK const &K, VectorN &fun) const
Definition RungeKutta.hh:1045
bool solve(VectorX const &t_mesh, VectorN const &ics, Solution< Real, N, M > &sol)
Definition RungeKutta.hh:1284
void max_projection_iterations(Integer t_max_projection_iterations)
Definition RungeKutta.hh:581
Eigen::Matrix< Real, N *S, N *S > MatrixJ
Definition RungeKutta.hh:55
VectorS b() const
Definition RungeKutta.hh:196
typename Tableau< Real, S >::Type Type
Definition RungeKutta.hh:72
void absolute_tolerance(Real t_absolute_tolerance)
Definition RungeKutta.hh:414
NewtonX m_newtonX
Definition RungeKutta.hh:76
bool is_dirk() const
Definition RungeKutta.hh:148
VectorS c() const
Definition RungeKutta.hh:208
void reverse(bool t_reverse)
Definition RungeKutta.hh:546
bool m_reverse
Definition RungeKutta.hh:89
bool m_adaptive
Definition RungeKutta.hh:87
Eigen::Vector< Real, Eigen::Dynamic > VectorX
Definition RungeKutta.hh:51
std::string name() const
Definition RungeKutta.hh:172
Real m_min_step
Definition RungeKutta.hh:85
void safety_factor(Real t_safety_factor)
Definition RungeKutta.hh:438
void implicit_system(typename ImplicitWrapper< Real, N, M >::FunctionF F, typename ImplicitWrapper< Real, N, M >::FunctionJF JF_x, typename ImplicitWrapper< Real, N, M >::FunctionJF JF_x_dot, typename ImplicitWrapper< Real, N, M >::FunctionH h=ImplicitWrapper< Real, N, M >::DefaultH, typename ImplicitWrapper< Real, N, M >::FunctionJH Jh_x=ImplicitWrapper< Real, N, M >::DefaultJH, typename ImplicitWrapper< Real, N, M >::FunctionID in_domain=ImplicitWrapper< Real, N, M >::DefaultID)
Definition RungeKutta.hh:231
void semi_explicit_system(std::string name, typename SemiExplicitWrapper< Real, N, M >::FunctionA A, typename SemiExplicitWrapper< Real, N, M >::FunctionTA TA_x, typename SemiExplicitWrapper< Real, N, M >::FunctionB b, typename SemiExplicitWrapper< Real, N, M >::FunctionJB Jb_x, typename SemiExplicitWrapper< Real, N, M >::FunctionH h=SemiExplicitWrapper< Real, N, M >::DefaultH, typename SemiExplicitWrapper< Real, N, M >::FunctionJH Jh_x=SemiExplicitWrapper< Real, N, M >::DefaultJH, typename SemiExplicitWrapper< Real, N, M >::FunctionID in_domain=SemiExplicitWrapper< Real, N, M >::DefaultID)
Definition RungeKutta.hh:384
void enable_adaptive_mode()
Definition RungeKutta.hh:503
Real projection_tolerance()
Definition RungeKutta.hh:562
bool project(VectorN const &x, Real t, VectorN &x_projected)
Definition RungeKutta.hh:1421
bool dirk_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hh:1119
void dirk_jacobian(Integer n, VectorN const &x, Real t, Real h, MatrixK const &K, MatrixN &jac) const
Definition RungeKutta.hh:1090
void verbose_mode(bool t_verbose)
Definition RungeKutta.hh:520
Real relative_tolerance()
Definition RungeKutta.hh:420
bool verbose_mode()
Definition RungeKutta.hh:514
void linear_system(typename LinearWrapper< Real, N, M >::FunctionE E, typename LinearWrapper< Real, N, M >::FunctionA A, typename LinearWrapper< Real, N, M >::FunctionB b, typename LinearWrapper< Real, N, M >::FunctionH h=LinearWrapper< Real, N, M >::DefaultH, typename LinearWrapper< Real, N, M >::FunctionJH Jh_x=LinearWrapper< Real, N, M >::DefaultJH, typename LinearWrapper< Real, N, M >::FunctionID in_domain=LinearWrapper< Real, N, M >::DefaultID)
Definition RungeKutta.hh:315
const Real SQRT_EPSILON
Definition RungeKutta.hh:69
Real estimate_order(std::vector< VectorX > const &t_mesh, VectorN const &ics, std::function< MatrixX(VectorX)> &sol)
Definition RungeKutta.hh:1560
void irk_jacobian(VectorN const &x, Real t, Real h, VectorK const &K, MatrixJ &jac) const
Definition RungeKutta.hh:938
void explicit_system(typename ExplicitWrapper< Real, N, M >::FunctionF f, typename ExplicitWrapper< Real, N, M >::FunctionJF Jf_x, typename ExplicitWrapper< Real, N, M >::FunctionH h=ExplicitWrapper< Real, N, M >::DefaultH, typename ExplicitWrapper< Real, N, M >::FunctionJH Jh_x=ExplicitWrapper< Real, N, M >::DefaultJH, typename ExplicitWrapper< Real, N, M >::FunctionID in_domain=ExplicitWrapper< Real, N, M >::DefaultID)
Definition RungeKutta.hh:274
Eigen::FullPivLU< MatrixP > m_lu
Definition RungeKutta.hh:91
void min_step(Real t_min_step)
Definition RungeKutta.hh:474
System m_system
Definition RungeKutta.hh:79
Tableau< Real, S > m_tableau
Definition RungeKutta.hh:78
Eigen::Matrix< Real, N, S > MatrixK
Definition RungeKutta.hh:54
Integer stages() const
Definition RungeKutta.hh:166
void max_substeps(Integer t_max_substeps)
Definition RungeKutta.hh:486
void enable_verbose_mode()
Definition RungeKutta.hh:529
Integer & max_projection_iterations()
Definition RungeKutta.hh:575
bool is_irk() const
Definition RungeKutta.hh:142
Optimist::RootFinder::Newton< Real, N > NewtonX
Definition RungeKutta.hh:58
bool erk_explicit_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new) const
Definition RungeKutta.hh:705
void projection_tolerance(Real t_projection_tolerance)
Definition RungeKutta.hh:568
Real m_projection_tolerance
Definition RungeKutta.hh:92
bool erk_implicit_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hh:822
bool advance(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hh:1194
void system(System t_system)
Definition RungeKutta.hh:220
void disable_adaptive_mode()
Definition RungeKutta.hh:508
Tableau< Real, S > const & tableau() const
Definition RungeKutta.hh:160
Real m_max_safety_factor
Definition RungeKutta.hh:84
void erk_implicit_function(Integer s, VectorN const &x, Real t, Real h, MatrixK const &K, VectorN &fun) const
Definition RungeKutta.hh:754
bool m_verbose
Definition RungeKutta.hh:88
Real m_safety_factor
Definition RungeKutta.hh:82
void disable_reverse_mode()
Definition RungeKutta.hh:556
Eigen::Vector< Real, N *S > VectorK
Definition RungeKutta.hh:53
void info(std::ostream &os)
Definition RungeKutta.hh:681
void relative_tolerance(Real t_relative_tolerance)
Definition RungeKutta.hh:426
void min_safety_factor(Real t_min_safety_factor)
Definition RungeKutta.hh:450
Type type() const
Definition RungeKutta.hh:130
void disable_projection()
Definition RungeKutta.hh:604
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Definition RungeKutta.hh:52
bool adaptive_solve(VectorX const &t_mesh, VectorN const &ics, Solution< Real, N, M > &sol)
Definition RungeKutta.hh:1346
Eigen::Vector< Real, Eigen::Dynamic > Time
Definition RungeKutta.hh:73
void disable_verbose_mode()
Definition RungeKutta.hh:534
Integer m_max_substeps
Definition RungeKutta.hh:86
bool is_embedded() const
Definition RungeKutta.hh:184
bool irk_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
Definition RungeKutta.hh:991
std::string info() const
Definition RungeKutta.hh:654
Real estimate_step(VectorN const &x, VectorN const &x_e, Real h_k) const
Definition RungeKutta.hh:640
void explicit_system(std::string name, typename ExplicitWrapper< Real, N, M >::FunctionF f, typename ExplicitWrapper< Real, N, M >::FunctionJF Jf_x, typename ExplicitWrapper< Real, N, M >::FunctionH h=ExplicitWrapper< Real, N, M >::DefaultH, typename ExplicitWrapper< Real, N, M >::FunctionJH Jh_x=ExplicitWrapper< Real, N, M >::DefaultJH, typename ExplicitWrapper< Real, N, M >::FunctionID in_domain=ExplicitWrapper< Real, N, M >::DefaultID)
Definition RungeKutta.hh:294
bool has_system()
Definition RungeKutta.hh:402
void erk_implicit_jacobian(Integer s, VectorN const &x, Real t, Real h, MatrixK const &K, MatrixN &jac) const
Definition RungeKutta.hh:797
typename Implicit< Real, N, M >::VectorH VectorM
Definition RungeKutta.hh:64
MatrixS A() const
Definition RungeKutta.hh:190
Real absolute_tolerance()
Definition RungeKutta.hh:408
Optimist::RootFinder::Newton< Real, N *S > NewtonK
Definition RungeKutta.hh:59
Integer m_max_projection_iterations
Definition RungeKutta.hh:93
NewtonK m_newtonK
Definition RungeKutta.hh:77
void implicit_system(std::string name, typename ImplicitWrapper< Real, N, M >::FunctionF F, typename ImplicitWrapper< Real, N, M >::FunctionJF JF_x, typename ImplicitWrapper< Real, N, M >::FunctionJF JF_x_dot, typename ImplicitWrapper< Real, N, M >::FunctionH h=ImplicitWrapper< Real, N, M >::DefaultH, typename ImplicitWrapper< Real, N, M >::FunctionJH Jh_x=ImplicitWrapper< Real, N, M >::DefaultJH, typename ImplicitWrapper< Real, N, M >::FunctionID in_domain=ImplicitWrapper< Real, N, M >::DefaultID)
Definition RungeKutta.hh:253
Real & max_safety_factor()
Definition RungeKutta.hh:456
bool m_projection
Definition RungeKutta.hh:94
std::function< VectorH(VectorF const &, Real)> FunctionH
Definition SemiExplicit.hh:278
static const FunctionH DefaultH
Definition SemiExplicit.hh:282
std::function< VectorB(VectorF const &, Real)> FunctionB
Definition SemiExplicit.hh:276
static const FunctionJH DefaultJH
Definition SemiExplicit.hh:283
std::function< MatrixA(VectorF const &, Real)> FunctionA
Definition SemiExplicit.hh:274
std::function< MatrixJB(VectorF const &, Real)> FunctionJB
Definition SemiExplicit.hh:277
static const FunctionID DefaultID
Definition SemiExplicit.hh:284
std::function< MatrixJH(VectorF const &, Real)> FunctionJH
Definition SemiExplicit.hh:279
std::function< bool(VectorF const &, Real)> FunctionID
Definition SemiExplicit.hh:280
std::function< TensorTA(VectorF const &, Real)> FunctionTA
Definition SemiExplicit.hh:275
The namespace for the Sandals library.
Definition Sandals.hh:89
SANDALS_DEFAULT_INTEGER_TYPE Integer
The Integer type as used for the API.
Definition Sandals.hh:97
Class container for the numerical solution of a system of ODEs/DAEs.
Definition Solution.hh:60
MatrixN x
Definition Solution.hh:66
Integer size() const
Definition Solution.hh:124
void resize(Integer size)
Definition Solution.hh:85
Vector t
Definition Solution.hh:65
void conservative_resize(Integer size)
Definition Solution.hh:95
MatrixM h
Definition Solution.hh:67
Struct container for the Butcher tableau of a Runge-Kutta method.
Definition Tableau.hh:38
enum class type :Integer {ERK=0, IRK=1, DIRK=2} Type
Definition Tableau.hh:42
Type type
Definition Tableau.hh:47
Integer order
Definition Tableau.hh:48
std::string name
Definition Tableau.hh:46
Vector b_e
Definition Tableau.hh:52
Eigen::Matrix< Real, S, S > Matrix
Definition Tableau.hh:44
Matrix A
Definition Tableau.hh:50
Eigen::Vector< Real, S > Vector
Definition Tableau.hh:43
Vector c
Definition Tableau.hh:53
Vector b
Definition Tableau.hh:51
bool is_embedded
Definition Tableau.hh:54