differential_equations/methods/erk/adaptive/
mod.rs

1//! Runge-Kutta solvers with support for dense output, embedded error estimation, and fixed steps.
2
3mod delay;
4mod ordinary;
5
6use crate::{
7    methods::{Adaptive, ExplicitRungeKutta},
8    tableau::ButcherTableau,
9    traits::{CallBackData, Real, State},
10};
11
12// Macro for adaptive step constructors
13macro_rules! impl_erk_adaptive_step_constructor {
14    ($method_name:ident, $fsal_val:expr, $order_val:expr, $s_val:expr, $m_val:expr, $doc:expr) => {
15        impl<E, T: Real, Y: State<T>, D: CallBackData>
16            ExplicitRungeKutta<E, Adaptive, T, Y, D, $order_val, $s_val, $m_val>
17        {
18            #[doc = $doc]
19            pub fn $method_name() -> Self {
20                let tableau = ButcherTableau::<T, $s_val, $m_val>::$method_name();
21                let c = tableau.c;
22                let a = tableau.a;
23                let b = tableau.b;
24                let bh = tableau.bh;
25                let bi = tableau.bi;
26                let fsal = $fsal_val;
27
28                ExplicitRungeKutta {
29                    c,
30                    a,
31                    b,
32                    bh,
33                    bi,
34                    fsal,
35                    ..Default::default()
36                }
37            }
38        }
39    };
40}
41
42// Adaptive step methods (embedded error estimation, cubic Hermite interpolation)
43impl_erk_adaptive_step_constructor!(
44    rkf45,
45    false,
46    5,
47    6,
48    6,
49    "Creates a Runge-Kutta-Fehlberg 4(5) method with error estimation."
50);
51impl_erk_adaptive_step_constructor!(
52    cash_karp,
53    false,
54    5,
55    6,
56    6,
57    "Creates a Cash-Karp 4(5) method with error estimation."
58);
59impl_erk_adaptive_step_constructor!(
60    rkv655e,
61    true,
62    6,
63    9,
64    10,
65    "Creates an Explicit Runge-Kutta 6(5) method with 9 stages and a 5th order interpolant."
66);
67impl_erk_adaptive_step_constructor!(
68    rkv656e,
69    true,
70    6,
71    9,
72    12,
73    "Creates an Explicit Runge-Kutta 6(5) method with 9 stages and a 6th order interpolant."
74);
75impl_erk_adaptive_step_constructor!(
76    rkv766e,
77    false,
78    7,
79    10,
80    13,
81    "Creates an Explicit Runge-Kutta 7(6) method with 10 stages and a 6th order interpolant."
82);
83impl_erk_adaptive_step_constructor!(
84    rkv767e,
85    false,
86    7,
87    10,
88    16,
89    "Creates an Explicit Runge-Kutta 7(6) method with 10 stages and a 7th order interpolant."
90);
91impl_erk_adaptive_step_constructor!(
92    rkv877e,
93    false,
94    8,
95    13,
96    17,
97    "Creates an Explicit Runge-Kutta 8(7) method with 13 stages and a 7th order interpolant."
98);
99impl_erk_adaptive_step_constructor!(
100    rkv878e,
101    false,
102    8,
103    13,
104    21,
105    "Creates an Explicit Runge-Kutta 8(7) method with 13 stages and an 8th order interpolant."
106);
107impl_erk_adaptive_step_constructor!(
108    rkv988e,
109    false,
110    9,
111    16,
112    21,
113    "Creates an Explicit Runge-Kutta 9(8) method with 16 stages and an 8th order interpolant."
114);
115impl_erk_adaptive_step_constructor!(
116    rkv989e,
117    false,
118    9,
119    16,
120    26,
121    "Creates an Explicit Runge-Kutta 9(8) method with 16 stages and a 9th order interpolant."
122);