Sandals  v0.0.0
A C++ library for ODEs/DAEs integration
Loading...
Searching...
No Matches
RungeKutta.hxx
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_HXX
14#define SANDALS_RUNGEKUTTA_HXX
15
16namespace Sandals {
17
18 /*\
19 | ____ _ __ _ _
20 | | _ \ _ _ _ __ __ _ ___| |/ / _| |_| |_ __ _
21 | | |_) | | | | '_ \ / _` |/ _ \ ' / | | | __| __/ _` |
22 | | _ <| |_| | | | | (_| | __/ . \ |_| | |_| || (_| |
23 | |_| \_\\__,_|_| |_|\__, |\___|_|\_\__,_|\__|\__\__,_|
24 | |___/
25 \*/
26
37 template <Integer S, Integer N, Integer M>
39 {
40 using VectorK = Eigen::Vector<Real, N*S>;
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>;
47 using VectorS = typename Tableau<S>::Vector;
48 using MatrixS = typename Tableau<S>::Matrix;
53
54 public:
56 using Type = typename Tableau<S>::Type;
57 using Time = Eigen::Vector<Real, Eigen::Dynamic>;
58
59 private:
71 bool m_adaptive{true};
72 bool m_verbose{false};
73 bool m_reverse{false};
74
75 Eigen::FullPivLU<MatrixP> m_lu;
78 bool m_projection{true};
79
80 public:
81
85 RungeKutta(const RungeKutta &) = delete;
86
90 RungeKutta & operator=(RungeKutta const &) = delete;
91
96 RungeKutta(Tableau<S> const &t_tableau)
97 : m_tableau(t_tableau) {
98 this->verbose_mode(this->m_verbose);
99 }
100
106 RungeKutta(Tableau<S> const &t_tableau, System t_system)
107 : m_tableau(t_tableau), m_system(t_system) {
108 this->verbose_mode(this->m_verbose);
109 }
110
115 Type type() const {return this->m_tableau.type;}
116
121 bool is_erk() const {return this->m_tableau.type == Type::ERK;}
122
127 bool is_irk() const {return this->m_tableau.type == Type::IRK;}
128
133 bool is_dirk() const {return this->m_tableau.type == Type::DIRK;}
134
139 Tableau<S> & tableau() {return this->m_tableau;}
140
145 Tableau<S> const & tableau() const {return this->m_tableau;}
146
151 Integer stages() const {return S;}
152
157 std::string name() const {return this->m_tableau.name;}
158
163 Integer order() const {return this->m_tableau.order;}
164
169 bool is_embedded() const {return this->m_tableau.is_embedded;}
170
175 MatrixS A() const {return this->m_tableau.A;}
176
181 VectorS b() const {return this->m_tableau.b;}
182
187 VectorS b_embedded() const {return this->m_tableau.b_e;}
188
193 VectorS c() const {return c;}
194
199 System system() {return this->m_system;}
200
205 void system(System t_system) {this->m_system = t_system;}
206
211 bool has_system() {return this->m_system != nullptr;}
212
218
223 void absolute_tolerance(Real t_absolute_tolerance) {this->m_absolute_tolerance = t_absolute_tolerance;}
224
230
235 void relative_tolerance(Real t_relative_tolerance) {this->m_relative_tolerance = t_relative_tolerance;}
236
242
247 void safety_factor(Real t_safety_factor) {this->m_safety_factor = t_safety_factor;}
248
254
259 void min_safety_factor(Real t_min_safety_factor) {this->m_min_safety_factor = t_min_safety_factor;}
260
266
271 void max_safety_factor(Real t_max_safety_factor) {this->m_max_safety_factor = t_max_safety_factor;}
272
277 Real & min_step() {return this->m_min_step;}
278
283 void min_step(Real t_min_step) {this->m_min_step = t_min_step;}
284
290
295 void max_substeps(Integer t_max_substeps) {this->m_max_substeps = t_max_substeps;}
296
301 bool adaptive_mode() {return this->m_adaptive;}
302
307 void adaptive(bool t_adaptive) {this->m_adaptive = t_adaptive;}
308
312 void enable_adaptive_mode() {this->m_adaptive = true;}
313
317 void disable_adaptive_mode() {this->m_adaptive = false;}
318
323 bool verbose_mode() {return this->m_verbose;}
324
329 void verbose_mode(bool t_verbose) {
330 this->m_verbose = t_verbose;
331 this->m_newtonX.verbose_mode(t_verbose);
332 this->m_newtonK.verbose_mode(t_verbose);
333 }
334
338 void enable_verbose_mode() {this->verbose_mode(true);}
339
343 void disable_verbose_mode() {this->verbose_mode(false);}
344
349 bool reverse_mode() {return this->m_reverse;}
350
355 void reverse(bool t_reverse) {this->m_reverse = t_reverse;}
356
360 void enable_reverse_mode() {this->m_reverse = true;}
361
365 void disable_reverse_mode() {this->m_reverse = false;}
366
372
377 void projection_tolerance(Real t_projection_tolerance)
378 {this->m_projection_tolerance = t_projection_tolerance;}
379
384 Integer & max_projection_iterations() {return this->m_max_projection_iterationsations;}
385
390 void max_projection_iterations(Integer t_max_projection_iterations)
391 {this->m_max_projection_iterationsations = t_max_projection_iterations;}
392
397 bool projection() {return this->m_projection;}
398
403 void projection(bool t_projection) {this->m_projection = t_projection;}
404
408 void enable_projection() {this->m_projection = true;}
409
413 void disable_projection() {this->m_projection = false;}
414
449 Real estimate_step(VectorN const &x, VectorN const &x_e, Real h_k) const
450 {
451 Real desired_error{this->m_absolute_tolerance + this->m_relative_tolerance *
452 std::max(x.array().abs().maxCoeff(), x_e.array().abs().maxCoeff())};
453 Real truncation_error{(x - x_e).array().abs().maxCoeff()};
454 return h_k * std::min(this->m_max_safety_factor, std::max(this->m_min_safety_factor,
455 this->m_safety_factor * std::pow(desired_error/truncation_error,
456 1.0/std::max(this->m_tableau.order, this->m_tableau.order_e))));
457 }
458
463 std::string info() const {
464 std::ostringstream os;
465 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
469 << "\t- type:\t";
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;
474 }
475 os
476 << std::endl
477 << "\t- embedded:\t" << this->is_embedded() << std::endl;
478 if (this->has_system()) {
479 os << "\t- system:\t" << this->m_system->name() << std::endl;
480 } else {
481 os << "\t- system:\t" << "none" << std::endl;
482 }
483 return os.str();
484 }
485
490 void info(std::ostream &os) {os << this->info();}
491
492 /*\
493 | _____ ____ _ __
494 | | ____| _ \| |/ /
495 | | _| | |_) | ' /
496 | | |___| _ <| . \
497 | |_____|_| \_\_|\_\
498 |
499 \*/
500
514 bool erk_explicit_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new) const
515 {
516 using Eigen::all;
517 using Eigen::seqN;
518
519 // Compute the K variables in the case of an explicit method and explicit system
520 VectorN x_node;
521 MatrixK K;
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();
524 if (!this->m_reverse) {
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));
526 } else {
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));
528 }
529 }
530 if (!K.allFinite()) {return false;}
531
532 // Perform the step and obtain the next state
533 x_new = x_old + K * this->m_tableau.b;
534
535 // Adapt next step
536 if (this->m_adaptive && this->m_tableau.is_embedded) {
537 VectorN x_emb = x_old + K * this->m_tableau.b_e;
538 h_new = this->estimate_step(x_new, x_emb, h_old);
539 }
540 return true;
541 }
542
563 void erk_implicit_function(Integer s, VectorN const &x, Real t, Real h, MatrixK const &K, VectorN &fun) const
564 {
565 using Eigen::all;
566 using Eigen::seqN;
567 VectorN x_node(x + K(all, seqN(0, s)) * this->m_tableau.A(s, seqN(0, s)).transpose());
568 if (!this->m_reverse) {
569 fun = this->m_system->F(x_node, K.col(s)/h, t + h * this->m_tableau.c(s));
570 } else {
571 fun = this->m_system->F_reverse(x_node, K.col(s)/h, t + h * this->m_tableau.c(s));
572 }
573 }
574
606 void erk_implicit_jacobian(Integer s, VectorN const &x, Real t, Real h, MatrixK const &K, MatrixN &jac) const
607 {
608 using Eigen::all;
609 using Eigen::seqN;
610 VectorN x_node(x + K(all, seqN(0, s)) * this->m_tableau.A(s, seqN(0, s)).transpose());
611 if (!this->m_reverse) {
612 jac = this->m_system->JF_x_dot(x_node, K.col(s)/h, t + h * this->m_tableau.c(s)) / h;
613 } else {
614 jac = this->m_system->JF_x_dot_reverse(x_node, K.col(s)/h, t + h * this->m_tableau.c(s)) / h;
615 }
616 }
617
631 bool erk_implicit_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
632 {
633 MatrixK K;
634 VectorN K_sol;
635 VectorN K_ini(VectorN::Zero());
636
637 // Check if the solver converged
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);},
644 K_ini, K_sol)) {
645 K.col(s) = K_sol;
646 } else {
647 return false;
648 }
649 }
650
651 // Perform the step and obtain the next state
652 x_new = x_old + K * this->m_tableau.b;
653
654 // Adapt next step
655 if (this->m_adaptive && this->m_tableau.is_embedded) {
656 VectorN x_emb(x_old + K * this->m_tableau.b_e);
657 h_new = this->estimate_step(x_new, x_emb, h_old);
658 }
659 return true;
660 }
661
662 /*\
663 | ___ ____ _ __
664 | |_ _| _ \| |/ /
665 | | || |_) | ' /
666 | | || _ <| . \
667 | |___|_| \_\_|\_\
668 |
669 \*/
670
695 void irk_function(VectorN const &x, Real t, Real h, VectorK const &K, VectorK &fun) const
696 {
697 VectorN x_node;
698 MatrixK K_mat{K.reshaped(N, S)};
699 MatrixK fun_mat;
700 for (Integer i{0}; i < S; ++i) {
701 x_node = x + K_mat * this->m_tableau.A.row(i).transpose();
702 if (!this->m_reverse) {
703 fun_mat.col(i) = this->m_system->F(x_node, K_mat.col(i)/h, t + h * this->m_tableau.c(i));
704 } else {
705 fun_mat.col(i) = this->m_system->F_reverse(x_node, K_mat.col(i)/h, t + h * this->m_tableau.c(i));
706 }
707 }
708 fun = fun_mat.reshaped(N*S, 1);
709 }
710
747 void irk_jacobian(VectorN const &x, Real t, Real h, VectorK const &K, MatrixJ &jac) const
748 {
749 using Eigen::seqN;
750
751 // Reset the Jacobian matrix
752 jac.setZero();
753
754 // Loop through each equation of the system
755 MatrixK K_mat{K.reshaped(N, S)};
756 Real t_node;
757 VectorN x_node, x_dot_node;
758 MatrixN JF_x, JF_x_dot;
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();
763
764 // Compute the Jacobians with respect to x and x_dot
765 x_dot_node = K_mat.col(i) / h;
766 if (!this->m_reverse) {
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);
769 } else {
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);
772 }
773
774 // Combine the Jacobians with respect to x and x_dot to obtain the Jacobian with respect to K
775 idx = seqN(i*N, N);
776 for (Integer j{0}; j < S; ++j) {
777 jdx = seqN(j*N, N);
778 if (i == j) {
779 jac(idx, jdx) = this->m_tableau.A(i,j) * JF_x + JF_x_dot / h;
780 } else {
781 jac(idx, jdx) = this->m_tableau.A(i,j) * JF_x;
782 }
783 }
784 }
785 }
786
800 bool irk_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
801 {
802 VectorK K;
803 VectorK K_ini(VectorK::Zero());
804
805 // Check if the solver converged
806 if (!this->m_newtonK.solve(
807 [this, &x_old, t_old, h_old](VectorK const &K_fun, VectorK &fun)
808 {this->irk_function(x_old, t_old, h_old, K_fun, fun);},
809 [this, &x_old, t_old, h_old](VectorK const &K_jac, MatrixJ &jac)
810 {this->irk_jacobian(x_old, t_old, h_old, K_jac, jac);},
811 K_ini, K))
812 {return false;}
813
814 // Perform the step and obtain the next state
815 x_new = x_old + K.reshaped(N, S) * this->m_tableau.b;
816
817 // Adapt next step
818 if (this->m_adaptive && this->m_tableau.is_embedded) {
819 VectorN x_emb(x_old + K.reshaped(N, S) * this->m_tableau.b_e);
820 h_new = this->estimate_step(x_new, x_emb, h_old);
821 }
822 return true;
823 }
824
825 /*\
826 | ____ ___ ____ _ __
827 | | _ \_ _| _ \| |/ /
828 | | | | | || |_) | ' /
829 | | |_| | || _ <| . \
830 | |____/___|_| \_\_|\_\
831 |
832 \*/
833
854 void dirk_function(Integer n, VectorN const &x, Real t, Real h, MatrixK const &K, VectorN &fun) const
855 {
856 using Eigen::all;
857 using Eigen::seqN;
858 VectorN x_node(x + K(all, seqN(0, n+1)) * this->m_tableau.A(n, seqN(0, n+1)).transpose());
859 if (!this->m_reverse) {
860 fun = this->m_system->F(x_node, K.col(n)/h, t + h * this->m_tableau.c(n));
861 } else {
862 fun = this->m_system->F_reverse(x_node, K.col(n)/h, t + h * this->m_tableau.c(n));
863 }
864 }
865
899 void dirk_jacobian(Integer n, VectorN const &x, Real t, Real h, MatrixK const &K, MatrixN &jac) const
900 {
901 using Eigen::all;
902 using Eigen::seqN;
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);
906 if (!this->m_reverse) {
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;
909 } else {
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;
912 }
913 }
914
928 bool dirk_step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
929 {
930 MatrixK K;
931 VectorN K_sol;
932 VectorN K_ini(VectorN::Zero());
933
934 // Check if the solver converged at each step
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);},
941 K_ini, K_sol)) {
942 K.col(n) = K_sol;
943 } else {
944 return false;
945 }
946 }
947
948 // Perform the step and obtain the next state
949 x_new = x_old + K * this->m_tableau.b;
950
951 // Adapt next step
952 if (this->m_adaptive && this->m_tableau.is_embedded) {
953 VectorN x_emb(x_old + K * this->m_tableau.b_e);
954 h_new = this->estimate_step(x_new, x_emb, h_old);
955 }
956 return true;
957 }
958
970 bool step(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
971 {
972 #define CMD "Sandals::RungeKutta::step(...): "
973
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.");
976
977 if (this->is_erk() && this->m_system->is_explicit()) {
978 return this->erk_explicit_step(x_old, t_old, h_old, x_new, h_new);
979 } else if (this->is_erk() && this->m_system->is_implicit()) {
980 return this->erk_implicit_step(x_old, t_old, h_old, x_new, h_new);
981 } else if (this->is_dirk()) {
982 return this->dirk_step(x_old, t_old, h_old, x_new, h_new);
983 } else {
984 return this->irk_step(x_old, t_old, h_old, x_new, h_new);
985 }
986
987 #undef CMD
988 }
989
1003 bool advance(VectorN const &x_old, Real t_old, Real h_old, VectorN &x_new, Real &h_new)
1004 {
1005 #define CMD "Sandals::RungeKutta::advance(...): "
1006
1007 // Check step size
1008 SANDALS_ASSERT(h_old > Real(0.0), CMD "in " << this->m_tableau.name << " solver, h = "<<
1009 h_old << ", expected > 0.");
1010
1011 // If the integration step failed, try again with substepping
1012 if (!this->step(x_old, t_old, h_old, x_new, h_new))
1013 {
1014 VectorN x_tmp(x_old);
1015 Real t_tmp{t_old}, h_tmp{h_old / Real(2.0)};
1016
1017 // Substepping logic
1018 Integer max_k{this->m_max_substeps * this->m_max_substeps}, k{2};
1019 Real h_new_tmp;
1020 while (k > 0) {
1021 // Calculate the next step with substepping logic
1022 if (this->step(x_tmp, t_tmp, h_tmp, x_new, h_new_tmp)) {
1023
1024 // Accept the step
1025 h_tmp = h_new_tmp;
1026
1027 // If substepping is enabled, double the step size
1028 if (k > 0 && k < max_k) {
1029 k -= 1;
1030 // If the substepping index is even, double the step size
1031 if (k % 2 == 0) {
1032 h_tmp = Real(2.0) * h_tmp;
1033 if (this->m_verbose) {
1034 SANDALS_WARNING(CMD "in " << this->m_tableau.name << " solver, at t = " << t_tmp <<
1035 ", integration succedded disable one substepping layer.");
1036 }
1037 }
1038 }
1039
1040 // Check the infinity norm of the solution
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.");
1043
1044 } else {
1045
1046 // If the substepping index is too high, abort the integration
1047 k += 2;
1048 SANDALS_ASSERT(k < max_k, CMD "in " << this->m_tableau.name << " solver, at t = " <<
1049 t_tmp << ", integration failed with h = " << h_tmp << ", aborting.");
1050 return false;
1051
1052 // Otherwise, try again with a smaller step
1053 if (this->m_verbose) {SANDALS_WARNING(CMD "in " << this->m_tableau.name << " solver, " <<
1054 "at t = " << t_tmp << ", integration failed, adding substepping layer.");}
1055 h_tmp /= Real(2.0);
1056 continue;
1057 }
1058
1059 // Store independent variable (or time)
1060 t_tmp += h_tmp;
1061 }
1062
1063 // Store output states substepping solutions
1064 x_new = x_tmp;
1065 h_new = h_tmp;
1066 }
1067
1068 // Project intermediate solution on the invariants
1069 if (this->m_projection) {
1070 VectorN x_projected;
1071 if (this->project(x_new, t_old + h_new, x_projected)) {
1072 x_new = x_projected;
1073 return true;
1074 } else {
1075 return false;
1076 }
1077 } else {
1078 return true;
1079 }
1080
1081 #undef CMD
1082 }
1083
1093 bool solve(VectorX const &t_mesh, VectorN const &ics, Solution<N, M> &sol)
1094 {
1095 using Eigen::last;
1096
1097 // Instantiate output
1098 sol.resize(t_mesh.size());
1099
1100 // Store first step
1101 sol.t(0) = t_mesh(0);
1102 sol.x.col(0) = ics;
1103 sol.h.col(0) = this->m_system->h(ics, t_mesh(0));
1104
1105 // Update the current step
1106 Integer step{0};
1107 VectorN x_step{ics};
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;
1110
1111 while (true) {
1112 // Integrate system
1113 if (!this->advance(sol.x.col(step), t_step, h_step, x_step, h_new_step)) {return false;}
1114
1115 // Update the current step
1116 t_step += h_step;
1117
1118 // Saturate the suggested timestep
1119 mesh_point_bool = std::abs(t_step - t_mesh(step+1)) < SQRT_EPSILON;
1120 saturation_bool = t_step + h_new_step > t_mesh(step+1) + SQRT_EPSILON;
1121 if (this->m_adaptive && this->m_tableau.is_embedded && !mesh_point_bool && saturation_bool) {
1122 h_tmp_step = h_new_step; // Used to store the previous step and keep the integration pace
1123 h_step = t_mesh(step+1) - t_step;
1124 } else {
1125 h_step = h_new_step;
1126 }
1127
1128 // Store solution if the step is a mesh point
1129 if (!this->m_adaptive || mesh_point_bool) {
1130 // Update temporaries
1131 step += 1;
1132 h_step = h_tmp_step;
1133
1134 // Update outputs
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);
1138
1139 // Check if the current step is the last one
1140 if (std::abs(t_step - t_mesh(last)) < SQRT_EPSILON) {break;}
1141 }
1142 }
1143 return true;
1144 }
1145
1155 bool adaptive_solve(VectorX const &t_mesh, VectorN const &ics, Solution<N, M> &sol)
1156 {
1157 using Eigen::all;
1158 using Eigen::last;
1159
1160 #define CMD "Sandals::RungeKutta::adaptive_solve(...): "
1161
1162 // Check if the adaptive method is enabled and the method is embedded
1163 if (!this->is_embedded()) {
1164 SANDALS_WARNING(CMD "the method is not embedded, using solve(...) method.");
1165 return this->solve(t_mesh, ics, sol);
1166 } else if (!this->m_adaptive) {
1167 SANDALS_WARNING(CMD "adaptive method is disabled, using solve(...) method.");
1168 return this->solve(t_mesh, ics, sol);
1169 }
1170
1171 // Instantiate output
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};
1174 if (this->m_tableau.is_embedded) {
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);
1177 } else {
1178 sol.resize(t_mesh.size());
1179 }
1180
1181 // Store first step
1182 sol.t(0) = t_mesh(0);
1183 sol.x.col(0) = ics;
1184 sol.h.col(0) = this->m_system->h(ics, t_mesh(0));
1185
1186 // Instantiate temporary variables
1187 Integer step{0};
1188 VectorN x_step{ics};
1189
1190 while (true) {
1191 // Integrate system
1192 this->advance(sol.x.col(step), sol.t(step), h_step, x_step, h_new_step);
1193
1194 // Saturate the suggested timestep
1195 if (this->m_adaptive && this->m_tableau.is_embedded) {
1196 h_step = std::max(std::min(h_new_step, h_max), h_min);
1197 }
1198
1199 SANDALS_ASSERT(step < sol.size(), CMD "safety length exceeded.");
1200
1201 // Store solution
1202 sol.t(step+1) = sol.t(step) + h_step;
1203 sol.x.col(step+1) = x_step;
1204 sol.h.col(step+1) = this->m_system->h(x_step, sol.t(step+1));
1205
1206 // Check if the current step is the last one
1207 if (sol.t(step+1) + h_step > t_mesh(last)) {break;}
1208
1209 // Update steps counter
1210 step += 1;
1211 }
1212
1213 // Resize the output
1215
1216 return true;
1217
1218 #undef CMD
1219 }
1220
1230 bool project(VectorN const &x, Real t, VectorN &x_projected)
1231 {
1232 #define CMD "Sandals::RungeKutta::project(...): "
1233
1234 // Check if there are any constraints
1235 x_projected = x;
1236 if (M > Integer(0)) {
1237
1238 VectorM h;
1239 MatrixM Jh_x;
1240 VectorP b, x_step;
1241 MatrixP A;
1242 A.setZero();
1243 A.template block<N, N>(0, 0) = MatrixN::Identity();
1244 for (Integer k{0}; k < this->m_max_projection_iterations; ++k) {
1245
1246 /* Standard projection method
1247 [A] {x} = {b}
1248 / I Jh_x^T \ / dx \ / x_t - x_k \
1249 | | | | = | |
1250 \ Jh_x 0 / \ lambda / \ -h /
1251 */
1252
1253 // Evaluate the invariants vector and its Jacobian
1254 h = this->m_system->h(x_projected, t);
1255 Jh_x = this->m_system->Jh_x(x_projected, t);
1256
1257 // Check if the solution is found
1258 if (h.norm() < this->m_projection_tolerance) {return true;}
1259
1260 // Build the linear system
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;
1265
1266 // Compute the solution of the linear system
1267 this->m_lu.compute(A);
1268 SANDALS_ASSERT(this->m_lu.rank() == N+M, CMD "singular Jacobian detected.");
1269 x_step = this->m_lu.solve(b);
1270
1271 // Check if the step is too small
1272 if (x_step.norm() < this->m_projection_tolerance * this->m_projection_tolerance) {return false;}
1273
1274 // Update the solution
1275 x_projected.noalias() += x_step(Eigen::seqN(0, N));
1276 }
1277 if (this->m_verbose) {SANDALS_WARNING(CMD "maximum number of iterations reached.");}
1278 return false;
1279 } else {
1280 return true;
1281 }
1282
1283 #undef CMD
1284 }
1285
1297 bool project_ics(VectorN const &x, Real t, std::vector<Integer> const & projected_equations,
1298 std::vector<Integer> const & projected_invariants, VectorN &x_projected) const
1299 {
1300 #define CMD "Sandals::RungeKutta::project_ics(...): "
1301
1302 Integer X{static_cast<Integer>(projected_equations.size())};
1303 Integer H{static_cast<Integer>(projected_invariants.size())};
1304
1305 // Check if there are any constraints
1306 x_projected = x;
1307 if (H > Integer(0)) {
1308
1309 VectorM h;
1310 MatrixM Jh_x;
1311 VectorX b(X+H), x_step(X+H);
1312 MatrixX A(X+H, X+H);
1313 A.setZero();
1314 A.block(0, 0, X, X) = MatrixX::Identity(X+H, X+H);
1315 Eigen::FullPivLU<MatrixX> lu;
1316 for (Integer k{0}; k < this->m_max_projection_iterations; ++k) {
1317
1318 /* Standard projection method
1319 [A] {x} = {b}
1320 / I Jh_x^T \ / dx \ / x_t - x_k \
1321 | | | | = | |
1322 \ Jh_x 0 / \ lambda / \ -h /
1323 */
1324
1325 // Evaluate the invariants vector and its Jacobian
1326 h = this->m_system->h(x_projected, t);
1327 Jh_x = this->m_system->Jh_x(x_projected, t);
1328
1329 // Select only the projected invariants
1330 h = h(projected_invariants);
1331 Jh_x = Jh_x(projected_invariants, projected_equations);
1332
1333 // Check if the solution is found
1334 if (h.norm() < this->m_projection_tolerance) {return true;}
1335
1336 // Build the linear system
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);
1340 b.tail(H) = -h;
1341
1342 // Compute the solution of the linear system
1343 lu.compute(A);
1344 SANDALS_ASSERT(lu.rank() == X+H, CMD "singular Jacobian detected.");
1345 x_step = this->m_lu.solve(b);
1346
1347 // Check if the step is too small
1348 if (x_step.norm() < this->m_projection_tolerance * this->m_projection_tolerance) {return false;}
1349
1350 // Update the solution
1351 x_projected(projected_equations).noalias() += x_step;
1352 }
1353 if (this->m_verbose) {SANDALS_WARNING(CMD "maximum number of iterations reached.");}
1354 return false;
1355 } else {
1356 return true;
1357 }
1358
1359 #undef CMD
1360 }
1361
1369 Real estimate_order(std::vector<VectorX> const &t_mesh, VectorN const &ics, std::function<MatrixX(VectorX)> &sol)
1370 {
1371 using Eigen::last;
1372
1373 #define CMD "Sandals::RungeKutta::estimate_order(...): "
1374
1375 SANDALS_ASSERT(t_mesh.size() > Integer(1), CMD "expected at least two time meshes.");
1376
1377 for (Integer i{0}; i < static_cast<Integer>(t_mesh.size()); ++i) {
1378
1379 // Check if the time mesh has the same initial and final time
1380 SANDALS_ASSERT((t_mesh[0](0) - t_mesh[i](0)) < SQRT_EPSILON,
1381 CMD "expected the same initial time.");
1382 SANDALS_ASSERT((t_mesh[0](last) - t_mesh[i](last)) < SQRT_EPSILON,
1383 CMD "expected the same final time.");
1384
1385 // Check if the mesh step is fixed all over the time mesh
1386 for (Integer j{1}; j < static_cast<Integer>(t_mesh[i].size()); ++j) {
1387 SANDALS_ASSERT((t_mesh[i](j) - t_mesh[i](j-1)) - (t_mesh[i](1) - t_mesh[i](0)) < SQRT_EPSILON,
1388 CMD "expected a fixed step.");
1389 }
1390 }
1391
1392 // Solve the system for each time scale
1393 Solution<N, M> sol_num;
1394 MatrixX sol_ana;
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) {
1397 SANDALS_ASSERT(this->solve(t_mesh[i], ics, sol_num), CMD "failed to solve the system for " <<
1398 "the" << i << "-th time mesh.");
1399 sol_ana = sol(sol_num.t);
1400 SANDALS_ASSERT(sol_ana.rows() == sol_num.x.rows(),
1401 CMD "expected the same number of states in analytical solution.");
1402 SANDALS_ASSERT(sol_ana.cols() == sol_num.x.cols(),
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();
1406 }
1407
1408 // Compute the order of the method thorugh least squares
1409 VectorX A(h_vec.array().log());
1410 VectorX b(e_vec.array().log());
1411 return ((A.transpose() * A).ldlt().solve(A.transpose() * b))(0);
1412
1413 #undef CMD
1414 }
1415
1416 }; // class RungeKutta
1417
1418} // namespace Sandals
1419
1420#endif // SANDALS_RUNGEKUTTA_HXX
#define CMD
#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