1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
//! Fixed Runge-Kutta methods for ODEs
use crate::{
error::Error,
interpolate::{Interpolation, cubic_hermite_interpolate},
methods::{ExplicitRungeKutta, Fixed, Ordinary},
ode::{ODE, OrdinaryNumericalMethod},
stats::Evals,
status::Status,
traits::{Real, State},
utils::validate_step_size_parameters,
};
impl<T: Real, Y: State<T>, const O: usize, const S: usize, const I: usize>
OrdinaryNumericalMethod<T, Y> for ExplicitRungeKutta<Ordinary, Fixed, T, Y, O, S, I>
{
fn init<F>(&mut self, ode: &F, t0: T, tf: T, y0: &Y) -> Result<Evals, Error<T, Y>>
where
F: ODE<T, Y>,
{
let mut evals = Evals::new();
// If h0 is zero, calculate initial step size for fixed-step methods
if self.h0 == T::zero() {
// Simple default step size for fixed-step methods
let duration = (tf - t0).abs();
let default_steps = T::from_usize(100).unwrap();
self.h0 = duration / default_steps;
}
// Check bounds
match validate_step_size_parameters::<T, Y>(self.h0, self.h_min, self.h_max, t0, tf) {
Ok(h0) => self.h = h0,
Err(status) => return Err(status),
} // Initialize Statistics
// Initialize State
self.t = t0;
self.y = *y0;
ode.diff(self.t, &self.y, &mut self.dydt);
evals.function += 1;
// Initialize previous state
self.t_prev = self.t;
self.y_prev = self.y;
self.dydt_prev = self.dydt;
// Initialize Status
self.status = Status::Initialized;
Ok(evals)
}
fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>
where
F: ODE<T, Y>,
{
let mut evals = Evals::new();
// Check max steps
if self.steps >= self.max_steps {
self.status = Status::Error(Error::MaxSteps {
t: self.t,
y: self.y,
});
return Err(Error::MaxSteps {
t: self.t,
y: self.y,
});
}
self.steps += 1;
// Save k[0] as the current derivative
self.k[0] = self.dydt;
// Compute stages
for i in 1..self.stages {
let mut y_stage = self.y;
for j in 0..i {
y_stage += self.k[j] * (self.a[i][j] * self.h);
}
ode.diff(self.t + self.c[i] * self.h, &y_stage, &mut self.k[i]);
}
evals.function += self.stages - 1; // We already have k[0]
// Store current state before update for interpolation
self.t_prev = self.t;
self.y_prev = self.y;
self.dydt_prev = self.k[0];
// Compute solution
let mut y_next = self.y;
for i in 0..self.stages {
y_next += self.k[i] * (self.b[i] * self.h);
}
// If method has dense output stages, compute them
if self.bi.is_some() {
// Compute extra stages for dense output
for i in 0..(I - S) {
let mut y_stage = self.y;
for j in 0..self.stages + i {
y_stage += self.k[j] * (self.a[self.stages + i][j] * self.h);
}
ode.diff(
self.t + self.c[self.stages + i] * self.h,
&y_stage,
&mut self.k[self.stages + i],
);
}
evals.function += I - S;
}
// Update state
self.t += self.h;
self.y = y_next;
// Calculate new derivative for next step
if self.fsal {
// If FSAL (First Same As Last) is enabled, we can reuse the last derivative
self.dydt = self.k[S - 1];
} else {
// Otherwise, compute the new derivative
ode.diff(self.t, &self.y, &mut self.dydt);
evals.function += 1;
}
self.status = Status::Solving;
Ok(evals)
}
fn t(&self) -> T {
self.t
}
fn y(&self) -> &Y {
&self.y
}
fn t_prev(&self) -> T {
self.t_prev
}
fn y_prev(&self) -> &Y {
&self.y_prev
}
fn h(&self) -> T {
self.h
}
fn set_h(&mut self, h: T) {
self.h = h;
}
fn status(&self) -> &Status<T, Y> {
&self.status
}
fn set_status(&mut self, status: Status<T, Y>) {
self.status = status;
}
}
impl<T: Real, Y: State<T>, const O: usize, const S: usize, const I: usize> Interpolation<T, Y>
for ExplicitRungeKutta<Ordinary, Fixed, T, Y, O, S, I>
{
fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>> {
// Check if t is within bounds
if t_interp < self.t_prev || t_interp > self.t {
return Err(Error::OutOfBounds {
t_interp,
t_prev: self.t_prev,
t_curr: self.t,
});
}
// If method has dense output coefficients, use them
if self.bi.is_some() {
// Calculate the normalized distance within the step [0, 1]
let s = (t_interp - self.t_prev) / self.h_prev;
// Get the interpolation coefficients
let bi = self.bi.as_ref().unwrap();
let mut cont = [T::zero(); I];
// Compute the interpolation coefficients using Horner's method
for i in 0..self.dense_stages {
// Start with the highest-order term
cont[i] = bi[i][self.order - 1];
// Apply Horner's method
for j in (0..self.order - 1).rev() {
cont[i] = cont[i] * s + bi[i][j];
}
// Multiply by s
cont[i] *= s;
}
// Compute the interpolated value
let mut y_interp = self.y_prev;
for i in 0..I {
y_interp += self.k[i] * cont[i] * self.h_prev;
}
Ok(y_interp)
} else {
// Otherwise use cubic Hermite interpolation
let y_interp = cubic_hermite_interpolate(
self.t_prev,
self.t,
&self.y_prev,
&self.y,
&self.dydt_prev,
&self.dydt,
t_interp,
);
Ok(y_interp)
}
}
}