conspire/math/integrate/dae/explicit/variable_step/verner_9/
mod.rs

1use crate::math::{
2    Scalar, Tensor, TensorVec,
3    integrate::{ExplicitDaeVariableStepExplicit, ode::explicit::variable_step::verner_9::*},
4};
5use std::ops::{Mul, Sub};
6
7impl<Y, Z, U, V> ExplicitDaeVariableStepExplicit<Y, Z, U, V> for Verner9
8where
9    Y: Tensor,
10    Z: Tensor,
11    U: TensorVec<Item = Y>,
12    V: TensorVec<Item = Z>,
13    for<'a> &'a Y: Mul<Scalar, Output = Y> + Sub<&'a Y, Output = Y>,
14{
15    fn slopes_solve(
16        mut evolution: impl FnMut(Scalar, &Y, &Z) -> Result<Y, String>,
17        mut solution: impl FnMut(Scalar, &Y, &Z) -> Result<Z, String>,
18        y: &Y,
19        z: &Z,
20        t: Scalar,
21        dt: Scalar,
22        k: &mut [Y],
23        y_trial: &mut Y,
24        z_trial: &mut Z,
25    ) -> Result<(), String> {
26        k[0] = evolution(t, y, z)?;
27        *y_trial = &k[0] * (A_2_1 * dt) + y;
28        *z_trial = solution(t + C_2 * dt, y_trial, z)?;
29        k[1] = evolution(t + C_2 * dt, y_trial, z_trial)?;
30        *y_trial = &k[0] * (A_3_1 * dt) + &k[1] * (A_3_2 * dt) + y;
31        *z_trial = solution(t + C_3 * dt, y_trial, z_trial)?;
32        k[2] = evolution(t + C_3 * dt, y_trial, z_trial)?;
33        *y_trial = &k[0] * (A_4_1 * dt) + &k[2] * (A_4_3 * dt) + y;
34        *z_trial = solution(t + C_4 * dt, y_trial, z_trial)?;
35        k[3] = evolution(t + C_4 * dt, y_trial, z_trial)?;
36        *y_trial = &k[0] * (A_5_1 * dt) + &k[2] * (A_5_3 * dt) + &k[3] * (A_5_4 * dt) + y;
37        *z_trial = solution(t + C_5 * dt, y_trial, z_trial)?;
38        k[4] = evolution(t + C_5 * dt, y_trial, z_trial)?;
39        *y_trial = &k[0] * (A_6_1 * dt) + &k[3] * (A_6_4 * dt) + &k[4] * (A_6_5 * dt) + y;
40        *z_trial = solution(t + C_6 * dt, y_trial, z_trial)?;
41        k[5] = evolution(t + C_6 * dt, y_trial, z_trial)?;
42        *y_trial = &k[0] * (A_7_1 * dt)
43            + &k[3] * (A_7_4 * dt)
44            + &k[4] * (A_7_5 * dt)
45            + &k[5] * (A_7_6 * dt)
46            + y;
47        *z_trial = solution(t + C_7 * dt, y_trial, z_trial)?;
48        k[6] = evolution(t + C_7 * dt, y_trial, z_trial)?;
49        *y_trial = &k[0] * (A_8_1 * dt) + &k[5] * (A_8_6 * dt) + &k[6] * (A_8_7 * dt) + y;
50        *z_trial = solution(t + C_8 * dt, y_trial, z_trial)?;
51        k[7] = evolution(t + C_8 * dt, y_trial, z_trial)?;
52        *y_trial = &k[0] * (A_9_1 * dt)
53            + &k[5] * (A_9_6 * dt)
54            + &k[6] * (A_9_7 * dt)
55            + &k[7] * (A_9_8 * dt)
56            + y;
57        *z_trial = solution(t + C_9 * dt, y_trial, z_trial)?;
58        k[8] = evolution(t + C_9 * dt, y_trial, z_trial)?;
59        *y_trial = &k[0] * (A_10_1 * dt)
60            + &k[5] * (A_10_6 * dt)
61            + &k[6] * (A_10_7 * dt)
62            + &k[7] * (A_10_8 * dt)
63            + &k[8] * (A_10_9 * dt)
64            + y;
65        *z_trial = solution(t + C_10 * dt, y_trial, z_trial)?;
66        k[9] = evolution(t + C_10 * dt, y_trial, z_trial)?;
67        *y_trial = &k[0] * (A_11_1 * dt)
68            + &k[5] * (A_11_6 * dt)
69            + &k[6] * (A_11_7 * dt)
70            + &k[7] * (A_11_8 * dt)
71            + &k[8] * (A_11_9 * dt)
72            + &k[9] * (A_11_10 * dt)
73            + y;
74        *z_trial = solution(t + C_11 * dt, y_trial, z_trial)?;
75        k[10] = evolution(t + C_11 * dt, y_trial, z_trial)?;
76        *y_trial = &k[0] * (A_12_1 * dt)
77            + &k[5] * (A_12_6 * dt)
78            + &k[6] * (A_12_7 * dt)
79            + &k[7] * (A_12_8 * dt)
80            + &k[8] * (A_12_9 * dt)
81            + &k[9] * (A_12_10 * dt)
82            + &k[10] * (A_12_11 * dt)
83            + y;
84        *z_trial = solution(t + C_12 * dt, y_trial, z_trial)?;
85        k[11] = evolution(t + C_12 * dt, y_trial, z_trial)?;
86        *y_trial = &k[0] * (A_13_1 * dt)
87            + &k[5] * (A_13_6 * dt)
88            + &k[6] * (A_13_7 * dt)
89            + &k[7] * (A_13_8 * dt)
90            + &k[8] * (A_13_9 * dt)
91            + &k[9] * (A_13_10 * dt)
92            + &k[10] * (A_13_11 * dt)
93            + &k[11] * (A_13_12 * dt)
94            + y;
95        *z_trial = solution(t + C_13 * dt, y_trial, z_trial)?;
96        k[12] = evolution(t + C_13 * dt, y_trial, z_trial)?;
97        *y_trial = &k[0] * (A_14_1 * dt)
98            + &k[5] * (A_14_6 * dt)
99            + &k[6] * (A_14_7 * dt)
100            + &k[7] * (A_14_8 * dt)
101            + &k[8] * (A_14_9 * dt)
102            + &k[9] * (A_14_10 * dt)
103            + &k[10] * (A_14_11 * dt)
104            + &k[11] * (A_14_12 * dt)
105            + &k[12] * (A_14_13 * dt)
106            + y;
107        *z_trial = solution(t + C_14 * dt, y_trial, z_trial)?;
108        k[13] = evolution(t + C_14 * dt, y_trial, z_trial)?;
109        *y_trial = &k[0] * (A_15_1 * dt)
110            + &k[5] * (A_15_6 * dt)
111            + &k[6] * (A_15_7 * dt)
112            + &k[7] * (A_15_8 * dt)
113            + &k[8] * (A_15_9 * dt)
114            + &k[9] * (A_15_10 * dt)
115            + &k[10] * (A_15_11 * dt)
116            + &k[11] * (A_15_12 * dt)
117            + &k[12] * (A_15_13 * dt)
118            + &k[13] * (A_15_14 * dt)
119            + y;
120        *z_trial = solution(t + dt, y_trial, z_trial)?;
121        k[14] = evolution(t + dt, y_trial, z_trial)?;
122        *y_trial = &k[0] * (A_16_1 * dt)
123            + &k[5] * (A_16_6 * dt)
124            + &k[6] * (A_16_7 * dt)
125            + &k[7] * (A_16_8 * dt)
126            + &k[8] * (A_16_9 * dt)
127            + &k[9] * (A_16_10 * dt)
128            + &k[10] * (A_16_11 * dt)
129            + &k[11] * (A_16_12 * dt)
130            + &k[12] * (A_16_13 * dt)
131            + y;
132        *z_trial = solution(t + dt, y_trial, z_trial)?;
133        k[15] = evolution(t + dt, y_trial, z_trial)?;
134        *y_trial = (&k[0] * B_1
135            + &k[7] * B_8
136            + &k[8] * B_9
137            + &k[9] * B_10
138            + &k[10] * B_11
139            + &k[11] * B_12
140            + &k[12] * B_13
141            + &k[13] * B_14
142            + &k[14] * B_15)
143            * dt
144            + y;
145        *z_trial = solution(t + dt, y_trial, z_trial)?;
146        Ok(())
147    }
148}