conspire/math/integrate/ode/explicit/variable_step/verner_8/
mod.rs1#[cfg(test)]
2mod test;
3
4use crate::math::{
5 Scalar, Tensor, TensorVec, Vector,
6 integrate::{Explicit, IntegrationError, OdeSolver, VariableStep, VariableStepExplicit},
7 interpolate::InterpolateSolution,
8};
9use crate::{ABS_TOL, REL_TOL};
10use std::ops::{Mul, Sub};
11
12const C_2: Scalar = 0.05;
13const C_3: Scalar = 0.1065625;
14const C_4: Scalar = 0.15984375;
15const C_5: Scalar = 0.39;
16const C_6: Scalar = 0.465;
17const C_7: Scalar = 0.155;
18const C_8: Scalar = 0.943;
19const C_9: Scalar = 0.901802041735857;
20const C_10: Scalar = 0.909;
21const C_11: Scalar = 0.94;
22
23const A_2_1: Scalar = 0.05;
24const A_3_1: Scalar = -0.0069931640625;
25const A_3_2: Scalar = 0.1135556640625;
26const A_4_1: Scalar = 0.0399609375;
27const A_4_3: Scalar = 0.1198828125;
28const A_5_1: Scalar = 0.36139756280045754;
29const A_5_3: Scalar = -1.3415240667004928;
30const A_5_4: Scalar = 1.3701265039000352;
31const A_6_1: Scalar = 0.049047202797202795;
32const A_6_4: Scalar = 0.23509720422144048;
33const A_6_5: Scalar = 0.18085559298135673;
34const A_7_1: Scalar = 0.06169289044289044;
35const A_7_4: Scalar = 0.11236568314640277;
36const A_7_5: Scalar = -0.03885046071451367;
37const A_7_6: Scalar = 0.01979188712522046;
38const A_8_1: Scalar = -1.767630240222327;
39const A_8_4: Scalar = -62.5;
40const A_8_5: Scalar = -6.061889377376669;
41const A_8_6: Scalar = 5.6508231982227635;
42const A_8_7: Scalar = 65.62169641937624;
43const A_9_1: Scalar = -1.1809450665549708;
44const A_9_4: Scalar = -41.50473441114321;
45const A_9_5: Scalar = -4.434438319103725;
46const A_9_6: Scalar = 4.260408188586133;
47const A_9_7: Scalar = 43.75364022446172;
48const A_9_8: Scalar = 0.00787142548991231;
49const A_10_1: Scalar = -1.2814059994414884;
50const A_10_4: Scalar = -45.047139960139866;
51const A_10_5: Scalar = -4.731362069449576;
52const A_10_6: Scalar = 4.514967016593808;
53const A_10_7: Scalar = 47.44909557172985;
54const A_10_8: Scalar = 0.01059228297111661;
55const A_10_9: Scalar = -0.0057468422638446166;
56const A_11_1: Scalar = -1.7244701342624853;
57const A_11_4: Scalar = -60.92349008483054;
58const A_11_5: Scalar = -5.951518376222392;
59const A_11_6: Scalar = 5.556523730698456;
60const A_11_7: Scalar = 63.98301198033305;
61const A_11_8: Scalar = 0.014642028250414961;
62const A_11_9: Scalar = 0.06460408772358203;
63const A_11_10: Scalar = -0.0793032316900888;
64const A_12_1: Scalar = -3.301622667747079;
65const A_12_4: Scalar = -118.01127235975251;
66const A_12_5: Scalar = -10.141422388456112;
67const A_12_6: Scalar = 9.139311332232058;
68const A_12_7: Scalar = 123.37594282840426;
69const A_12_8: Scalar = 4.62324437887458;
70const A_12_9: Scalar = -3.3832777380682018;
71const A_12_10: Scalar = 4.527592100324618;
72const A_12_11: Scalar = -5.828495485811623;
73const A_13_1: Scalar = -3.039515033766309;
74const A_13_4: Scalar = -109.26086808941763;
75const A_13_5: Scalar = -9.290642497400293;
76const A_13_6: Scalar = 8.43050498176491;
77const A_13_7: Scalar = 114.20100103783314;
78const A_13_8: Scalar = -0.9637271342145479;
79const A_13_9: Scalar = -5.0348840888021895;
80const A_13_10: Scalar = 5.958130824002923;
81
82const B_1: Scalar = 0.04427989419007951;
83const B_6: Scalar = 0.3541049391724449;
84const B_7: Scalar = 0.24796921549564377;
85const B_8: Scalar = -15.694202038838085;
86const B_9: Scalar = 25.084064965558564;
87const B_10: Scalar = -31.738367786260277;
88const B_11: Scalar = 22.938283273988784;
89const B_12: Scalar = -0.2361324633071542;
90
91const D_1: Scalar = -0.00003272103901028138;
92const D_6: Scalar = -0.0005046250618777704;
93const D_7: Scalar = 0.0001211723589784759;
94const D_8: Scalar = -20.142336771313868;
95const D_9: Scalar = 5.2371785994398286;
96const D_10: Scalar = -8.156744408794658;
97const D_11: Scalar = 22.938283273988784;
98const D_12: Scalar = -0.2361324633071542;
99const D_13: Scalar = 0.36016794372897754;
100
101#[doc = include_str!("doc.md")]
102#[derive(Debug)]
103pub struct Verner8 {
104 pub abs_tol: Scalar,
106 pub rel_tol: Scalar,
108 pub dt_beta: Scalar,
110 pub dt_expn: Scalar,
112 pub dt_cut: Scalar,
114 pub dt_min: Scalar,
116}
117
118impl Default for Verner8 {
119 fn default() -> Self {
120 Self {
121 abs_tol: ABS_TOL,
122 rel_tol: REL_TOL,
123 dt_beta: 0.9,
124 dt_expn: 8.0,
125 dt_cut: 0.5,
126 dt_min: ABS_TOL,
127 }
128 }
129}
130
131impl<Y, U> OdeSolver<Y, U> for Verner8
132where
133 Y: Tensor,
134 U: TensorVec<Item = Y>,
135{
136}
137
138impl VariableStep for Verner8 {
139 fn abs_tol(&self) -> Scalar {
140 self.abs_tol
141 }
142 fn rel_tol(&self) -> Scalar {
143 self.rel_tol
144 }
145 fn dt_beta(&self) -> Scalar {
146 self.dt_beta
147 }
148 fn dt_expn(&self) -> Scalar {
149 self.dt_expn
150 }
151 fn dt_cut(&self) -> Scalar {
152 self.dt_cut
153 }
154 fn dt_min(&self) -> Scalar {
155 self.dt_min
156 }
157}
158
159impl<Y, U> Explicit<Y, U> for Verner8
160where
161 Self: OdeSolver<Y, U>,
162 Y: Tensor,
163 for<'a> &'a Y: Mul<Scalar, Output = Y> + Sub<&'a Y, Output = Y>,
164 U: TensorVec<Item = Y>,
165{
166 const SLOPES: usize = 13;
167 fn integrate(
168 &self,
169 function: impl FnMut(Scalar, &Y) -> Result<Y, String>,
170 time: &[Scalar],
171 initial_condition: Y,
172 ) -> Result<(Vector, U, U), IntegrationError> {
173 self.integrate_variable_step(function, time, initial_condition)
174 }
175}
176
177impl<Y, U> VariableStepExplicit<Y, U> for Verner8
178where
179 Self: OdeSolver<Y, U>,
180 Y: Tensor,
181 for<'a> &'a Y: Mul<Scalar, Output = Y> + Sub<&'a Y, Output = Y>,
182 U: TensorVec<Item = Y>,
183{
184 fn slopes(
185 &self,
186 mut function: impl FnMut(Scalar, &Y) -> Result<Y, String>,
187 y: &Y,
188 t: Scalar,
189 dt: Scalar,
190 k: &mut [Y],
191 y_trial: &mut Y,
192 ) -> Result<Scalar, String> {
193 k[0] = function(t, y)?;
194 *y_trial = &k[0] * (A_2_1 * dt) + y;
195 k[1] = function(t + C_2 * dt, y_trial)?;
196 *y_trial = &k[0] * (A_3_1 * dt) + &k[1] * (A_3_2 * dt) + y;
197 k[2] = function(t + C_3 * dt, y_trial)?;
198 *y_trial = &k[0] * (A_4_1 * dt) + &k[2] * (A_4_3 * dt) + y;
199 k[3] = function(t + C_4 * dt, y_trial)?;
200 *y_trial = &k[0] * (A_5_1 * dt) + &k[2] * (A_5_3 * dt) + &k[3] * (A_5_4 * dt) + y;
201 k[4] = function(t + C_5 * dt, y_trial)?;
202 *y_trial = &k[0] * (A_6_1 * dt) + &k[3] * (A_6_4 * dt) + &k[4] * (A_6_5 * dt) + y;
203 k[5] = function(t + C_6 * dt, y_trial)?;
204 *y_trial = &k[0] * (A_7_1 * dt)
205 + &k[3] * (A_7_4 * dt)
206 + &k[4] * (A_7_5 * dt)
207 + &k[5] * (A_7_6 * dt)
208 + y;
209 k[6] = function(t + C_7 * dt, y_trial)?;
210 *y_trial = &k[0] * (A_8_1 * dt)
211 + &k[3] * (A_8_4 * dt)
212 + &k[4] * (A_8_5 * dt)
213 + &k[5] * (A_8_6 * dt)
214 + &k[6] * (A_8_7 * dt)
215 + y;
216 k[7] = function(t + C_8 * dt, y_trial)?;
217 *y_trial = &k[0] * (A_9_1 * dt)
218 + &k[3] * (A_9_4 * dt)
219 + &k[4] * (A_9_5 * dt)
220 + &k[5] * (A_9_6 * dt)
221 + &k[6] * (A_9_7 * dt)
222 + &k[7] * (A_9_8 * dt)
223 + y;
224 k[8] = function(t + C_9 * dt, y_trial)?;
225 *y_trial = &k[0] * (A_10_1 * dt)
226 + &k[3] * (A_10_4 * dt)
227 + &k[4] * (A_10_5 * dt)
228 + &k[5] * (A_10_6 * dt)
229 + &k[6] * (A_10_7 * dt)
230 + &k[7] * (A_10_8 * dt)
231 + &k[8] * (A_10_9 * dt)
232 + y;
233 k[9] = function(t + C_10 * dt, y_trial)?;
234 *y_trial = &k[0] * (A_11_1 * dt)
235 + &k[3] * (A_11_4 * dt)
236 + &k[4] * (A_11_5 * dt)
237 + &k[5] * (A_11_6 * dt)
238 + &k[6] * (A_11_7 * dt)
239 + &k[7] * (A_11_8 * dt)
240 + &k[8] * (A_11_9 * dt)
241 + &k[9] * (A_11_10 * dt)
242 + y;
243 k[10] = function(t + C_11 * dt, y_trial)?;
244 *y_trial = &k[0] * (A_12_1 * dt)
245 + &k[3] * (A_12_4 * dt)
246 + &k[4] * (A_12_5 * dt)
247 + &k[5] * (A_12_6 * dt)
248 + &k[6] * (A_12_7 * dt)
249 + &k[7] * (A_12_8 * dt)
250 + &k[8] * (A_12_9 * dt)
251 + &k[9] * (A_12_10 * dt)
252 + &k[10] * (A_12_11 * dt)
253 + y;
254 k[11] = function(t + dt, y_trial)?;
255 *y_trial = &k[0] * (A_13_1 * dt)
256 + &k[3] * (A_13_4 * dt)
257 + &k[4] * (A_13_5 * dt)
258 + &k[5] * (A_13_6 * dt)
259 + &k[6] * (A_13_7 * dt)
260 + &k[7] * (A_13_8 * dt)
261 + &k[8] * (A_13_9 * dt)
262 + &k[9] * (A_13_10 * dt)
263 + y;
264 k[12] = function(t + dt, y_trial)?;
265 *y_trial = (&k[0] * B_1
266 + &k[5] * B_6
267 + &k[6] * B_7
268 + &k[7] * B_8
269 + &k[8] * B_9
270 + &k[9] * B_10
271 + &k[10] * B_11
272 + &k[11] * B_12)
273 * dt
274 + y;
275 Ok(((&k[0] * D_1
276 + &k[5] * D_6
277 + &k[6] * D_7
278 + &k[7] * D_8
279 + &k[8] * D_9
280 + &k[9] * D_10
281 + &k[10] * D_11
282 + &k[11] * D_12
283 + &k[12] * D_13)
284 * dt)
285 .norm_inf())
286 }
287 fn step(
288 &self,
289 mut function: impl FnMut(Scalar, &Y) -> Result<Y, String>,
290 y: &mut Y,
291 t: &mut Scalar,
292 y_sol: &mut U,
293 t_sol: &mut Vector,
294 dydt_sol: &mut U,
295 dt: &mut Scalar,
296 _k: &mut [Y],
297 y_trial: &Y,
298 e: Scalar,
299 ) -> Result<(), String> {
300 if e < self.abs_tol || e / y_trial.norm_inf() < self.rel_tol {
301 *t += *dt;
302 *y = y_trial.clone();
303 t_sol.push(*t);
304 y_sol.push(y.clone());
305 dydt_sol.push(function(*t, y)?);
306 }
307 self.time_step(e, dt);
308 Ok(())
309 }
310}
311
312impl<Y, U> InterpolateSolution<Y, U> for Verner8
313where
314 Y: Tensor,
315 for<'a> &'a Y: Mul<Scalar, Output = Y> + Sub<&'a Y, Output = Y>,
316 U: TensorVec<Item = Y>,
317{
318 fn interpolate(
319 &self,
320 time: &Vector,
321 tp: &Vector,
322 yp: &U,
323 mut function: impl FnMut(Scalar, &Y) -> Result<Y, String>,
324 ) -> Result<(U, U), IntegrationError> {
325 let mut dt;
326 let mut i;
327 let mut k_1;
328 let mut k_2;
329 let mut k_3;
330 let mut k_4;
331 let mut k_5;
332 let mut k_6;
333 let mut k_7;
334 let mut k_8;
335 let mut k_9;
336 let mut k_10;
337 let mut k_11;
338 let mut k_12;
339 let mut t;
340 let mut y;
341 let mut y_int = U::new();
342 let mut dydt_int = U::new();
343 let mut y_trial;
344 for time_k in time.iter() {
345 i = tp.iter().position(|tp_i| tp_i >= time_k).unwrap();
346 if time_k == &tp[i] {
347 t = tp[i];
348 y_trial = yp[i].clone();
349 dt = 0.0;
350 } else {
351 t = tp[i - 1];
352 y = yp[i - 1].clone();
353 dt = time_k - t;
354 k_1 = function(t, &y)?;
355 y_trial = &k_1 * (A_2_1 * dt) + &y;
356 k_2 = function(t + C_2 * dt, &y_trial)?;
357 y_trial = &k_1 * (A_3_1 * dt) + &k_2 * (A_3_2 * dt) + &y;
358 k_3 = function(t + C_3 * dt, &y_trial)?;
359 y_trial = &k_1 * (A_4_1 * dt) + &k_3 * (A_4_3 * dt) + &y;
360 k_4 = function(t + C_4 * dt, &y_trial)?;
361 y_trial = &k_1 * (A_5_1 * dt) + &k_3 * (A_5_3 * dt) + &k_4 * (A_5_4 * dt) + &y;
362 k_5 = function(t + C_5 * dt, &y_trial)?;
363 y_trial = &k_1 * (A_6_1 * dt) + &k_4 * (A_6_4 * dt) + &k_5 * (A_6_5 * dt) + &y;
364 k_6 = function(t + C_6 * dt, &y_trial)?;
365 y_trial = &k_1 * (A_7_1 * dt)
366 + &k_4 * (A_7_4 * dt)
367 + &k_5 * (A_7_5 * dt)
368 + &k_6 * (A_7_6 * dt)
369 + &y;
370 k_7 = function(t + C_7 * dt, &y_trial)?;
371 y_trial = &k_1 * (A_8_1 * dt)
372 + &k_4 * (A_8_4 * dt)
373 + &k_5 * (A_8_5 * dt)
374 + &k_6 * (A_8_6 * dt)
375 + &k_7 * (A_8_7 * dt)
376 + &y;
377 k_8 = function(t + C_8 * dt, &y_trial)?;
378 y_trial = &k_1 * (A_9_1 * dt)
379 + &k_4 * (A_9_4 * dt)
380 + &k_5 * (A_9_5 * dt)
381 + &k_6 * (A_9_6 * dt)
382 + &k_7 * (A_9_7 * dt)
383 + &k_8 * (A_9_8 * dt)
384 + &y;
385 k_9 = function(t + C_9 * dt, &y_trial)?;
386 y_trial = &k_1 * (A_10_1 * dt)
387 + &k_4 * (A_10_4 * dt)
388 + &k_5 * (A_10_5 * dt)
389 + &k_6 * (A_10_6 * dt)
390 + &k_7 * (A_10_7 * dt)
391 + &k_8 * (A_10_8 * dt)
392 + &k_9 * (A_10_9 * dt)
393 + &y;
394 k_10 = function(t + C_10 * dt, &y_trial)?;
395 y_trial = &k_1 * (A_11_1 * dt)
396 + &k_4 * (A_11_4 * dt)
397 + &k_5 * (A_11_5 * dt)
398 + &k_6 * (A_11_6 * dt)
399 + &k_7 * (A_11_7 * dt)
400 + &k_8 * (A_11_8 * dt)
401 + &k_9 * (A_11_9 * dt)
402 + &k_10 * (A_11_10 * dt)
403 + &y;
404 k_11 = function(t + C_11 * dt, &y_trial)?;
405 y_trial = &k_1 * (A_12_1 * dt)
406 + &k_4 * (A_12_4 * dt)
407 + &k_5 * (A_12_5 * dt)
408 + &k_6 * (A_12_6 * dt)
409 + &k_7 * (A_12_7 * dt)
410 + &k_8 * (A_12_8 * dt)
411 + &k_9 * (A_12_9 * dt)
412 + &k_10 * (A_12_10 * dt)
413 + &k_11 * (A_12_11 * dt)
414 + &y;
415 k_12 = function(t + dt, &y_trial)?;
416 y_trial = (&k_1 * B_1
417 + &k_6 * B_6
418 + &k_7 * B_7
419 + &k_8 * B_8
420 + &k_9 * B_9
421 + &k_10 * B_10
422 + &k_11 * B_11
423 + &k_12 * B_12)
424 * dt
425 + &y;
426 }
427 dydt_int.push(function(t + dt, &y_trial)?);
428 y_int.push(y_trial);
429 }
430 Ok((y_int, dydt_int))
431 }
432}