Skip to main content

numra_ode/
sensitivity.rs

1//! Forward sensitivity analysis for parameterised ODE systems.
2//!
3//! Given an ODE system parameterised by a vector $p \in \mathbb{R}^{N_s}$,
4//!
5//! ```text
6//! dy/dt = f(t, y, p),    y(t_0) = y_0(p),
7//! ```
8//!
9//! the *forward sensitivity matrix* `S(t) = ∂y(t) / ∂p ∈ ℝ^{N × N_s}` satisfies
10//! the variational equations
11//!
12//! ```text
13//! dS/dt = J_y(t, y, p) · S + J_p(t, y, p),    S(t_0) = ∂y_0/∂p
14//! ```
15//!
16//! where `J_y = ∂f/∂y ∈ ℝ^{N×N}` and `J_p = ∂f/∂p ∈ ℝ^{N×N_s}`.
17//!
18//! Numra integrates the augmented state `z = [y; vec(S)] ∈ ℝ^{N(1+N_s)}` with
19//! the user's chosen [`crate::Solver`] (DoPri5, Tsit5, Vern\*, Radau5, BDF,
20//! Esdirk\*). The augmented Jacobian is `block_diag(J_y, …, J_y)` —
21//! the standard CVODES *simultaneous-corrector* formulation.
22//!
23//! # Layout conventions
24//!
25//! Numra mixes row-major and column-major flattening so that each Jacobian's
26//! "natural" axis is contiguous. Implementors of [`ParametricOdeSystem`] must
27//! follow these rules exactly:
28//!
29//! | Quantity                        | Layout         | Index map                         |
30//! | ------------------------------- | -------------- | --------------------------------- |
31//! | state `y`                       | flat, length `N`           | `y[i]`                |
32//! | state Jacobian `J_y`            | row-major, `N × N`         | `jy[i*N + j] = ∂f_i/∂y_j` |
33//! | parameter Jacobian `J_p`        | column-major, `N × N_s`    | `jp[k*N + i] = ∂f_i/∂p_k` |
34//! | sensitivity `S`                 | column-major, `N × N_s`    | `s[k*N + i] = ∂y_i/∂p_k`  |
35//! | augmented state `z`             | `[y; vec(S)]`              | `z[N + k*N + i] = S_{i,k}` |
36//! | initial sensitivity `∂y_0/∂p`   | column-major, `N × N_s`    | `s0[k*N + i] = ∂y_{0,i}/∂p_k` |
37//!
38//! Choosing column-major for `J_p` and `S` matches CVODES indexing, makes the
39//! per-parameter sub-vector `S_{:,k}` contiguous (and free to slice via
40//! [`SensitivityResult::sensitivity_for_param`]), and improves cache locality
41//! in the inner loop `Ṡ_{:,k} = J_y · S_{:,k} + J_p_{:,k}`. `J_y` stays
42//! row-major because state-row access is the standard `OdeSystem::jacobian`
43//! convention used by every implicit solver in this crate.
44//!
45//! # Worked example: 2-state, 2-parameter linear ODE
46//!
47//! Consider the system
48//!
49//! ```text
50//! ẏ_0 = -p_0 · y_0 + p_1 · y_1
51//! ẏ_1 =          - p_1 · y_1
52//! ```
53//!
54//! Then
55//!
56//! ```text
57//! J_y = [ -p_0    p_1 ]      J_p = [ -y_0    y_1 ]
58//!       [   0    -p_1 ]            [   0    -y_1 ]
59//! ```
60//!
61//! With `N = N_s = 2`, `J_y` is row-major and `J_p` is column-major.
62//! Implementing the trait by hand:
63//!
64//! ```
65//! use numra_ode::sensitivity::ParametricOdeSystem;
66//!
67//! struct Lin2 { p: [f64; 2] }
68//!
69//! impl ParametricOdeSystem<f64> for Lin2 {
70//!     fn n_states(&self) -> usize { 2 }
71//!     fn n_params(&self) -> usize { 2 }
72//!     fn params(&self) -> &[f64] { &self.p }
73//!
74//!     fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
75//!         dy[0] = -p[0] * y[0] + p[1] * y[1];
76//!         dy[1] =                -p[1] * y[1];
77//!     }
78//!
79//!     fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) {
80//!         // Row-major (N × N): jy[i*N + j] = ∂f_i/∂y_j
81//!         jy[0*2 + 0] = -self.p[0];   // ∂f_0/∂y_0
82//!         jy[0*2 + 1] =  self.p[1];   // ∂f_0/∂y_1
83//!         jy[1*2 + 0] =  0.0;         // ∂f_1/∂y_0
84//!         jy[1*2 + 1] = -self.p[1];   // ∂f_1/∂y_1
85//!     }
86//!
87//!     fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) {
88//!         // Column-major (N × N_s): jp[k*N + i] = ∂f_i/∂p_k
89//!         jp[0*2 + 0] = -y[0];        // ∂f_0/∂p_0
90//!         jp[0*2 + 1] =  0.0;         // ∂f_1/∂p_0
91//!         jp[1*2 + 0] =  y[1];        // ∂f_0/∂p_1
92//!         jp[1*2 + 1] = -y[1];        // ∂f_1/∂p_1
93//!     }
94//!
95//!     // Flag analytical so AugmentedSystem honors the overrides above
96//!     // instead of falling through to its inline-FD fast path.
97//!     fn has_analytical_jacobian_y(&self) -> bool { true }
98//!     fn has_analytical_jacobian_p(&self) -> bool { true }
99//! }
100//!
101//! let sys = Lin2 { p: [1.0, 0.5] };
102//! assert_eq!(sys.n_states(), 2);
103//! assert_eq!(sys.n_params(), 2);
104//! ```
105//!
106//! Note the asymmetry: `J_y[0*2 + 1]` is row 0, column 1 of `J_y` (state
107//! Jacobian, row-major), while `J_p[1*2 + 0]` is row 0, column 1 of `J_p`
108//! (parameter Jacobian, column-major). They look indexed the same way in
109//! Rust syntax but address different mathematical entries by design.
110//!
111//! # Send / Sync
112//!
113//! [`ParametricOdeSystem`] does **not** require `Send + Sync`, matching the
114//! [`crate::OdeSystem`] precedent. Add the bounds at the call site
115//! (`where Sys: ParametricOdeSystem<S> + Send + Sync`) when sharing the
116//! system across threads — for instance in parallel Monte Carlo. Tightening
117//! a trait bound is a breaking change; loosening it is not, so the default
118//! is permissive.
119//!
120//! # Performance note (v1)
121//!
122//! Numra's v1 ships *correct* forward sensitivity through every solver, with
123//! analytical `J_y` and `J_p` overrides replacing FD on the dominant cost
124//! path. The block-diagonal structure of the augmented Jacobian — which
125//! would let implicit solvers reuse a single LU of `M = (1/γh)·I_N - J_y`
126//! across all `N_s + 1` sub-systems — is **not yet exploited**: Radau5 and
127//! BDF currently perform dense LU on the full `N(N_s+1) × N(N_s+1)` matrix.
128//! Block-aware factorisation is tracked as a named follow-up (see
129//! `docs/internal-followups.md`).
130//!
131//! Author: Moussa Leblouba
132//! Date: 4 February 2026
133//! Modified: 6 May 2026
134
135use numra_core::Scalar;
136
137use crate::error::SolverError;
138use crate::problem::{OdeProblem, OdeSystem};
139use crate::solver::{Solver, SolverOptions, SolverStats};
140
141/// An ODE system parameterised by a parameter vector `p`.
142///
143/// Implementors expose dimensions, the nominal parameter vector, the
144/// parameterised right-hand side, and optionally analytical Jacobians.
145/// Forward finite differences are provided as defaults for `jacobian_y` and
146/// `jacobian_p`, so a minimal implementation supplies only `n_states`,
147/// `n_params`, `params`, and `rhs_with_params`.
148///
149/// See the [module-level documentation](self) for the layout conventions
150/// (`J_y` row-major, `J_p` and `S` column-major) and a worked 2-state,
151/// 2-parameter example.
152pub trait ParametricOdeSystem<S: Scalar> {
153    /// Number of state variables `N` (the dimension of `y`).
154    fn n_states(&self) -> usize;
155
156    /// Number of parameters `N_s` (the dimension of `p`).
157    fn n_params(&self) -> usize;
158
159    /// Nominal parameter vector. Length must equal [`Self::n_params`].
160    fn params(&self) -> &[S];
161
162    /// Evaluate `f(t, y, p) -> dy/dt` for an arbitrary parameter vector `p`.
163    ///
164    /// This is the primitive form: it lets the FD-default Jacobians perturb
165    /// `p` without mutating `self`.
166    fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]);
167
168    /// Convenience: evaluate the RHS at the system's nominal parameters.
169    fn rhs(&self, t: S, y: &[S], dydt: &mut [S]) {
170        self.rhs_with_params(t, y, self.params(), dydt);
171    }
172
173    /// Fill the state Jacobian `J_y[i,j] = ∂f_i/∂y_j` in **row-major** order
174    /// (`jac[i*N + j]`, length `N²`).
175    ///
176    /// Default: forward finite differences via [`Self::rhs_with_params`].
177    /// Override to supply an analytical Jacobian — required for performance
178    /// on stiff problems with implicit solvers.
179    fn jacobian_y(&self, t: S, y: &[S], jac: &mut [S]) {
180        let n = self.n_states();
181        let p = self.params();
182        let h_factor = S::EPSILON.sqrt();
183
184        let mut f0 = vec![S::ZERO; n];
185        let mut f1 = vec![S::ZERO; n];
186        let mut y_pert = y.to_vec();
187
188        self.rhs_with_params(t, y, p, &mut f0);
189
190        for j in 0..n {
191            let yj = y_pert[j];
192            let h = h_factor * (S::ONE + yj.abs());
193            y_pert[j] = yj + h;
194            self.rhs_with_params(t, &y_pert, p, &mut f1);
195            y_pert[j] = yj;
196            for i in 0..n {
197                jac[i * n + j] = (f1[i] - f0[i]) / h;
198            }
199        }
200    }
201
202    /// Fill the parameter Jacobian `J_p[i,k] = ∂f_i/∂p_k` in **column-major**
203    /// order (`jp[k*N + i]`, length `N · N_s`).
204    ///
205    /// Default: forward finite differences. Override for performance.
206    fn jacobian_p(&self, t: S, y: &[S], jp: &mut [S]) {
207        let n = self.n_states();
208        let np = self.n_params();
209        let p_nominal = self.params();
210        let h_factor = S::EPSILON.sqrt();
211
212        let mut f0 = vec![S::ZERO; n];
213        let mut f1 = vec![S::ZERO; n];
214        let mut p_pert = p_nominal.to_vec();
215
216        self.rhs_with_params(t, y, p_nominal, &mut f0);
217
218        for k in 0..np {
219            let pk = p_pert[k];
220            let h = h_factor * (S::ONE + pk.abs());
221            p_pert[k] = pk + h;
222            self.rhs_with_params(t, y, &p_pert, &mut f1);
223            p_pert[k] = pk;
224            for i in 0..n {
225                jp[k * n + i] = (f1[i] - f0[i]) / h;
226            }
227        }
228    }
229
230    /// Initial sensitivity `∂y_0/∂p`, column-major (`s0[k*N + i]`, length
231    /// `N · N_s`).
232    ///
233    /// Defaults to zero, which is correct when `y_0` does not depend on `p`.
234    /// Override only if `y_0 = y_0(p)`.
235    fn initial_sensitivity(&self, _y0: &[S], s0: &mut [S]) {
236        for s in s0.iter_mut() {
237            *s = S::ZERO;
238        }
239    }
240
241    /// Returns `true` iff [`Self::jacobian_y`] has been overridden with an
242    /// analytical implementation. Default: `false`.
243    ///
244    /// `AugmentedSystem` checks this flag inside its hot path. When `false`
245    /// (the default), it inlines its own forward-FD using reused scratch
246    /// buffers — bypassing both the trait default and any override.
247    /// When `true`, it calls `system.jacobian_y(...)` directly so the user's
248    /// analytical override is honored.
249    ///
250    /// This is the pattern used by CVODES (`set_jacobian_user_supplied`):
251    /// flagging analytical Jacobians lets the augmented-system path skip
252    /// FD scratch allocation in the common FD case while still respecting
253    /// analytical overrides in the stiff-problem case where they materially
254    /// matter. **If you override `jacobian_y`, also override this method to
255    /// return `true`** — otherwise your override is silently bypassed when
256    /// running through `solve_forward_sensitivity`.
257    ///
258    /// # Debug-build safety net
259    ///
260    /// To catch the silent-misconfiguration case (override the method,
261    /// forget the flag), `AugmentedSystem` runs a one-time consistency
262    /// check on the first RHS call: it evaluates [`Self::jacobian_y`]
263    /// (which dispatches to the user's override or the FD default) and
264    /// compares against the inline FD result. If they disagree by more
265    /// than a Frobenius-norm relative threshold of `1e-3`, the system
266    /// panics in debug builds with a message naming the missing flag
267    /// override. Release builds skip the check entirely (zero overhead).
268    /// This is a safety net, not a contract — fixing the panic by
269    /// returning `true` here is the canonical resolution.
270    fn has_analytical_jacobian_y(&self) -> bool {
271        false
272    }
273
274    /// Returns `true` iff [`Self::jacobian_p`] has been overridden with an
275    /// analytical implementation. Default: `false`. See
276    /// [`Self::has_analytical_jacobian_y`] for the rationale, contract,
277    /// and debug-build consistency check.
278    fn has_analytical_jacobian_p(&self) -> bool {
279        false
280    }
281
282    // ------------------------------------------------------------------
283    // Mass-matrix surface (F-IC-SENS-MASS-SURFACE, 0.1.5, Foundation Spec
284    // §6 #14 / §3.8). Mirrors [`OdeSystem`]'s mass-matrix surface so that
285    // parameter-sensitivity hosts can declare DAE structure that the
286    // augmented `OdeSystem`-impl on [`AugmentedSystem`] then lifts onto
287    // the augmented `(N · (1 + N_s))`-dimensional state.
288    //
289    // Defaults match [`OdeSystem`]'s defaults exactly: `is_autonomous =
290    // false`, `has_mass_matrix = false`, `mass_matrix` fills row-major
291    // identity, `is_singular_mass = false`, `algebraic_indices = empty`.
292    // Every pre-0.1.5 implementor compiles unchanged.
293    // ------------------------------------------------------------------
294
295    /// Is the system autonomous? (f does not depend on t explicitly)
296    ///
297    /// Default: `false`. Override to `true` when `f(t, y, p) = f(y, p)`
298    /// independent of `t`. `AugmentedSystem` forwards this from the host
299    /// so the augmented system inherits the host's time-dependence
300    /// structure — variational equations inherit autonomy from the
301    /// state equation.
302    fn is_autonomous(&self) -> bool {
303        false
304    }
305
306    /// Does this system have a mass matrix?
307    ///
308    /// Default: `false` (the standard parameterised ODE `dy/dt = f(t, y, p)`).
309    /// Override to `true` for DAE-shaped parameterised systems
310    /// `M · y' = f(t, y, p)` where `M` is non-identity or singular.
311    /// See [`OdeSystem::has_mass_matrix`] for the equivalent surface on
312    /// the unparameterised side.
313    fn has_mass_matrix(&self) -> bool {
314        false
315    }
316
317    /// Get the mass matrix `M` for the DAE: `M · y' = f(t, y, p)`.
318    ///
319    /// Default returns identity (standard ODE), regardless of what
320    /// [`Self::has_mass_matrix`] reports — callers should consult
321    /// [`Self::has_mass_matrix`] before relying on this output.
322    ///
323    /// Layout: row-major, `mass[i * n + j] = M[i, j]`, length `n²` where
324    /// `n = self.n_states()`. Mirrors [`OdeSystem::mass_matrix`].
325    fn mass_matrix(&self, mass: &mut [S]) {
326        let n = self.n_states();
327        for i in 0..n {
328            for j in 0..n {
329                mass[i * n + j] = if i == j { S::ONE } else { S::ZERO };
330            }
331        }
332    }
333
334    /// Is the mass matrix singular? (i.e., is this a DAE?)
335    ///
336    /// Default: `false`. Override to `true` for semi-explicit DAEs where
337    /// one or more rows of `M` are zero (algebraic constraints).
338    fn is_singular_mass(&self) -> bool {
339        false
340    }
341
342    /// Indices of algebraic variables (rows `i` where `M[i, i] = 0`).
343    ///
344    /// Default: empty (all rows differential). Override for DAE systems
345    /// to enumerate the algebraic rows so downstream solvers can apply
346    /// row-aware error scaling and the augmented-system lift can map
347    /// each host algebraic index `i` to its `(N_s + 1)` copies on the
348    /// augmented diagonal.
349    fn algebraic_indices(&self) -> Vec<usize> {
350        Vec::new()
351    }
352}
353
354/// Lift the host's `n × n` mass matrix `M` block-diagonally onto the
355/// augmented system's `(n + n·m) × (n + n·m)` state, where `m =
356/// n_params` is the number of sensitivity columns.
357///
358/// # Variational-equation derivation
359///
360/// The host DAE is `M · y' = f(t, y, p)`. Differentiating each side with
361/// respect to a parameter `p_k` gives the variational equation
362///
363/// ```text
364///   M · (∂y/∂p_k)' = J_y · (∂y/∂p_k) + J_p_{:,k}
365/// ```
366///
367/// so each sensitivity column `∂y/∂p_k` lives in the same M-weighted
368/// space as the state itself. Stacking the augmented state `z = [y;
369/// vec(S)] ∈ ℝ^{n(1+m)}` (column-major sensitivity flattening, see
370/// module-level layout conventions) and writing the augmented dynamics in
371/// matrix form `M_aug · z' = F(t, z, p)` therefore requires
372///
373/// ```text
374///   M_aug = block_diag(M, M, …, M)  with (m + 1) copies of M
375/// ```
376///
377/// — one M for the state block, plus one M for each of the `m`
378/// sensitivity-column blocks. Singular rows of `M` lift through to every
379/// block (see `lift_algebraic_indices_for_augmented`).
380///
381/// # Layout
382///
383/// `host_m`: row-major, length `n²`, `host_m[i * n + j] = M[i, j]`.
384/// `aug_mass`: row-major, length `(n * (1 + m))²`, filled by this
385/// function. All cross-block entries are zero (block-diagonal).
386fn lift_mass_matrix_for_augmented<S: Scalar>(host_m: &[S], n: usize, m: usize, aug_mass: &mut [S]) {
387    let dim = n * (1 + m);
388    debug_assert_eq!(host_m.len(), n * n);
389    debug_assert_eq!(aug_mass.len(), dim * dim);
390    // Zero the augmented matrix first (all off-block entries are zero).
391    for slot in aug_mass.iter_mut() {
392        *slot = S::ZERO;
393    }
394    // Place `(m + 1)` copies of M on the diagonal blocks. Block `b`
395    // occupies rows `[b*n, (b+1)*n)` and columns `[b*n, (b+1)*n)`.
396    for b in 0..=m {
397        let row0 = b * n;
398        let col0 = b * n;
399        for i in 0..n {
400            for j in 0..n {
401                aug_mass[(row0 + i) * dim + (col0 + j)] = host_m[i * n + j];
402            }
403        }
404    }
405}
406
407/// Lift the host's algebraic indices onto the augmented system's `(N_s
408/// + 1)` blocks.
409///
410/// The host's algebraic structure (rows `i` of `M` that are entirely
411/// zero) propagates to every sensitivity column through the variational
412/// equation: the row-`i` variational equation
413/// `M[i,:] · (∂y/∂p_k)' = J_y[i,:] · (∂y/∂p_k) + J_p[i,k]` is algebraic
414/// in every block whenever it is algebraic in the host. So for each host
415/// algebraic index `i`, the augmented algebraic indices are
416/// `{i, n + i, 2n + i, …, m · n + i}` — one copy per block.
417///
418/// Result is sorted ascending to match the conventional `algebraic_indices`
419/// output ordering.
420fn lift_algebraic_indices_for_augmented(host_indices: &[usize], n: usize, m: usize) -> Vec<usize> {
421    let mut out = Vec::with_capacity(host_indices.len() * (1 + m));
422    for b in 0..=m {
423        let offset = b * n;
424        for &i in host_indices {
425            out.push(offset + i);
426        }
427    }
428    out.sort_unstable();
429    out
430}
431
432/// Wraps a [`ParametricOdeSystem`] as an [`OdeSystem`] over the augmented
433/// state `z = [y; vec(S)]` (column-major sensitivity flattening).
434///
435/// The augmented dimension is `N · (1 + N_s)`. The augmented Jacobian is
436/// `block_diag(J_y, J_y, …, J_y)` (CVODES simultaneous-corrector form).
437pub struct AugmentedSystem<S: Scalar, Sys: ParametricOdeSystem<S>> {
438    /// The underlying parameterised system.
439    pub system: Sys,
440    jy_scratch: std::cell::RefCell<Vec<S>>,
441    jp_scratch: std::cell::RefCell<Vec<S>>,
442    // Reusable FD scratch lifted out of the trait's per-call allocations.
443    // Touched only on the FD path (when has_analytical_jacobian_{y,p}
444    // returns false). Each RefCell is borrowed in isolation; no overlap.
445    fd_f0: std::cell::RefCell<Vec<S>>,
446    fd_f1: std::cell::RefCell<Vec<S>>,
447    fd_y_pert: std::cell::RefCell<Vec<S>>,
448    fd_p_pert: std::cell::RefCell<Vec<S>>,
449    // Lifted host algebraic indices for the augmented `(N · (1 + N_s))`
450    // state. Computed once at construction by `lift_algebraic_indices_for_
451    // augmented` from `system.algebraic_indices()`. Algebraic indices are
452    // assumed constant over the host's lifetime (the trait contract does
453    // not require this explicitly, but every shipping implementor returns
454    // a constant set, and a future foundation change that ever makes them
455    // dynamic should revisit this cache).
456    augmented_algebraic_indices: Vec<usize>,
457    // Debug-only: tracks whether the first-call analytical-vs-FD
458    // consistency check has run. See `check_jacobian_flags`.
459    #[cfg(debug_assertions)]
460    flag_check_done: std::cell::Cell<bool>,
461}
462
463impl<S: Scalar, Sys: ParametricOdeSystem<S>> AugmentedSystem<S, Sys> {
464    /// Wrap a [`ParametricOdeSystem`] into an `OdeSystem`-shaped augmented
465    /// system suitable for any [`crate::Solver`].
466    pub fn new(system: Sys) -> Self {
467        let n = system.n_states();
468        let np = system.n_params();
469        // Snapshot the host's algebraic indices once at construction and
470        // lift them onto the augmented state's `(N_s + 1)` blocks. See
471        // `lift_algebraic_indices_for_augmented` for the lift derivation
472        // and the host-stability assumption.
473        let augmented_algebraic_indices =
474            lift_algebraic_indices_for_augmented(&system.algebraic_indices(), n, np);
475        Self {
476            system,
477            jy_scratch: std::cell::RefCell::new(vec![S::ZERO; n * n]),
478            jp_scratch: std::cell::RefCell::new(vec![S::ZERO; n * np]),
479            fd_f0: std::cell::RefCell::new(vec![S::ZERO; n]),
480            fd_f1: std::cell::RefCell::new(vec![S::ZERO; n]),
481            fd_y_pert: std::cell::RefCell::new(vec![S::ZERO; n]),
482            fd_p_pert: std::cell::RefCell::new(vec![S::ZERO; np]),
483            augmented_algebraic_indices,
484            #[cfg(debug_assertions)]
485            flag_check_done: std::cell::Cell::new(false),
486        }
487    }
488
489    /// Debug-build consistency check for the analytical-Jacobian flags.
490    ///
491    /// Catches the silent-misconfiguration class of bug where a user
492    /// implements `jacobian_y` (or `jacobian_p`) analytically but forgets
493    /// to override `has_analytical_jacobian_y` (or `_p`) — in which case
494    /// the augmented hot path bypasses the analytical override entirely
495    /// and runs FD instead. The user sees "Numra is slow on my problem"
496    /// with no diagnostic path.
497    ///
498    /// On the first `rhs` call (and only then), if a flag returns `false`,
499    /// we evaluate `system.jacobian_*` (which dispatches to the user's
500    /// override or the trait's FD default) and compare against the inline
501    /// FD result. If the user did not override the method, both calls
502    /// compute FD and agree to floating-point noise. If the user did
503    /// override with analytical, the two will typically disagree well
504    /// beyond the relative threshold — and we panic with a clear message.
505    ///
506    /// The threshold is `1e-3` relative (Frobenius). FD is only accurate
507    /// to `~sqrt(eps_mach) ≈ 1.5e-8`, so 1e-3 is loose enough to absorb
508    /// honest FD-vs-analytical disagreement at the initial state without
509    /// false positives, while still catching the "forgot the flag"
510    /// scenario where the analytical and FD results differ by O(1).
511    ///
512    /// Wrapped in `#[cfg(debug_assertions)]` — release builds pay zero.
513    #[cfg(debug_assertions)]
514    fn check_jacobian_flags(&self, t: S, y: &[S]) {
515        let n = self.system.n_states();
516        let np = self.system.n_params();
517        let threshold = S::from_f64(1e-3);
518
519        if !self.system.has_analytical_jacobian_y() {
520            let mut user_jy = vec![S::ZERO; n * n];
521            let mut fd_jy = vec![S::ZERO; n * n];
522            self.system.jacobian_y(t, y, &mut user_jy);
523            self.fd_jacobian_y_inline(t, y, &mut fd_jy);
524            if jacobians_differ_significantly(&user_jy, &fd_jy, threshold) {
525                panic!(
526                    "ParametricOdeSystem implements `jacobian_y` analytically \
527                     but `has_analytical_jacobian_y()` returns `false`; the \
528                     analytical implementation is being ignored and `AugmentedSystem` \
529                     is running finite differences instead. Override \
530                     `has_analytical_jacobian_y` to return `true`. \
531                     (This check is debug-build only; release builds will \
532                     silently use FD.)"
533                );
534            }
535        }
536
537        if !self.system.has_analytical_jacobian_p() && np > 0 {
538            let mut user_jp = vec![S::ZERO; n * np];
539            let mut fd_jp = vec![S::ZERO; n * np];
540            self.system.jacobian_p(t, y, &mut user_jp);
541            self.fd_jacobian_p_inline(t, y, &mut fd_jp);
542            if jacobians_differ_significantly(&user_jp, &fd_jp, threshold) {
543                panic!(
544                    "ParametricOdeSystem implements `jacobian_p` analytically \
545                     but `has_analytical_jacobian_p()` returns `false`; the \
546                     analytical implementation is being ignored and `AugmentedSystem` \
547                     is running finite differences instead. Override \
548                     `has_analytical_jacobian_p` to return `true`. \
549                     (This check is debug-build only; release builds will \
550                     silently use FD.)"
551                );
552            }
553        }
554    }
555
556    /// Forward-FD `J_y` using AugmentedSystem-local scratch buffers. Only
557    /// invoked when `system.has_analytical_jacobian_y()` is `false`.
558    fn fd_jacobian_y_inline(&self, t: S, y: &[S], jy: &mut [S]) {
559        let n = self.system.n_states();
560        let p = self.system.params();
561        let h_factor = S::EPSILON.sqrt();
562        let mut f0 = self.fd_f0.borrow_mut();
563        let mut f1 = self.fd_f1.borrow_mut();
564        let mut y_pert = self.fd_y_pert.borrow_mut();
565
566        self.system.rhs_with_params(t, y, p, &mut f0);
567        y_pert.copy_from_slice(y);
568        for j in 0..n {
569            let yj = y_pert[j];
570            let h = h_factor * (S::ONE + yj.abs());
571            y_pert[j] = yj + h;
572            self.system.rhs_with_params(t, &y_pert, p, &mut f1);
573            y_pert[j] = yj;
574            for i in 0..n {
575                jy[i * n + j] = (f1[i] - f0[i]) / h;
576            }
577        }
578    }
579
580    /// Forward-FD `J_p` using AugmentedSystem-local scratch buffers. Only
581    /// invoked when `system.has_analytical_jacobian_p()` is `false`.
582    fn fd_jacobian_p_inline(&self, t: S, y: &[S], jp: &mut [S]) {
583        let n = self.system.n_states();
584        let np = self.system.n_params();
585        let p_nominal = self.system.params();
586        let h_factor = S::EPSILON.sqrt();
587        let mut f0 = self.fd_f0.borrow_mut();
588        let mut f1 = self.fd_f1.borrow_mut();
589        let mut p_pert = self.fd_p_pert.borrow_mut();
590
591        self.system.rhs_with_params(t, y, p_nominal, &mut f0);
592        p_pert.copy_from_slice(p_nominal);
593        for k in 0..np {
594            let pk = p_pert[k];
595            let h = h_factor * (S::ONE + pk.abs());
596            p_pert[k] = pk + h;
597            self.system.rhs_with_params(t, y, &p_pert, &mut f1);
598            p_pert[k] = pk;
599            for i in 0..n {
600                jp[k * n + i] = (f1[i] - f0[i]) / h;
601            }
602        }
603    }
604
605    /// Total augmented dimension `N · (1 + N_s)`.
606    pub fn augmented_dim(&self) -> usize {
607        self.system.n_states() * (1 + self.system.n_params())
608    }
609
610    /// Build the initial augmented state `[y_0; vec(∂y_0/∂p)]`.
611    pub fn initial_augmented(&self, y0: &[S]) -> Vec<S> {
612        let n = self.system.n_states();
613        let np = self.system.n_params();
614
615        let mut z0 = Vec::with_capacity(n * (1 + np));
616        z0.extend_from_slice(y0);
617
618        let mut s0 = vec![S::ZERO; n * np];
619        self.system.initial_sensitivity(y0, &mut s0);
620        z0.extend_from_slice(&s0);
621        z0
622    }
623}
624
625/// Frobenius-norm relative comparison: returns `true` if
626/// `||a - b||_F / max(||b||_F, 1) > threshold`. Used by the debug-build
627/// analytical-vs-FD Jacobian flag check.
628#[cfg(debug_assertions)]
629fn jacobians_differ_significantly<S: Scalar>(a: &[S], b: &[S], threshold: S) -> bool {
630    debug_assert_eq!(a.len(), b.len());
631    let mut diff_sq = S::ZERO;
632    let mut b_sq = S::ZERO;
633    for i in 0..a.len() {
634        let d = a[i] - b[i];
635        diff_sq = diff_sq + d * d;
636        b_sq = b_sq + b[i] * b[i];
637    }
638    let denom = b_sq.sqrt().max(S::ONE);
639    let rel = diff_sq.sqrt() / denom;
640    rel > threshold
641}
642
643impl<S: Scalar, Sys: ParametricOdeSystem<S>> OdeSystem<S> for AugmentedSystem<S, Sys> {
644    fn dim(&self) -> usize {
645        self.augmented_dim()
646    }
647
648    /// Augmented RHS: `dy/dt = f(t,y,p)`, `dS_{:,k}/dt = J_y · S_{:,k} + J_p_{:,k}`.
649    fn rhs(&self, t: S, z: &[S], dz: &mut [S]) {
650        let n = self.system.n_states();
651        let np = self.system.n_params();
652        let y = &z[..n];
653
654        // Debug-only one-time consistency check: analytical-Jacobian flags
655        // vs the user's actual `jacobian_*` override. See
656        // `check_jacobian_flags` for rationale.
657        #[cfg(debug_assertions)]
658        if !self.flag_check_done.get() {
659            self.flag_check_done.set(true);
660            self.check_jacobian_flags(t, y);
661        }
662
663        // (a) original dynamics.
664        self.system.rhs(t, y, &mut dz[..n]);
665
666        // (b) Jacobians at the current state. Analytical when the system
667        // flags it; otherwise inline FD using local scratch buffers (lifted
668        // out of the trait defaults so the integration's hot path doesn't
669        // pay per-call allocation).
670        let mut jy = self.jy_scratch.borrow_mut();
671        let mut jp = self.jp_scratch.borrow_mut();
672        if self.system.has_analytical_jacobian_y() {
673            self.system.jacobian_y(t, y, &mut jy);
674        } else {
675            drop(jy);
676            self.fd_jacobian_y_inline(t, y, &mut self.jy_scratch.borrow_mut());
677            jy = self.jy_scratch.borrow_mut();
678        }
679        if self.system.has_analytical_jacobian_p() {
680            self.system.jacobian_p(t, y, &mut jp);
681        } else {
682            drop(jp);
683            self.fd_jacobian_p_inline(t, y, &mut self.jp_scratch.borrow_mut());
684            jp = self.jp_scratch.borrow_mut();
685        }
686
687        // (c) per-parameter sensitivity column: dS_{:,k}/dt = J_y · S_{:,k} + J_p_{:,k}.
688        // Sensitivity column-major: S_{i,k} lives at z[n + k*n + i].
689        for k in 0..np {
690            for i in 0..n {
691                let mut acc = S::ZERO;
692                // J_y row-major dot S_{:,k} (contiguous slice of z).
693                for j in 0..n {
694                    acc = acc + jy[i * n + j] * z[n + k * n + j];
695                }
696                acc = acc + jp[k * n + i];
697                dz[n + k * n + i] = acc;
698            }
699        }
700    }
701
702    /// Augmented Jacobian, row-major over the full `N(1+N_s) × N(1+N_s)` block.
703    ///
704    /// Fills the `N_s + 1` diagonal blocks with `J_y`. Off-diagonal blocks
705    /// (which would carry `∂(J_y · S_{:,k})/∂y` second-derivative coupling)
706    /// are zero — the standard CVODES *simultaneous-corrector* approximation.
707    /// Implicit solvers in v1 perform dense LU on the full block; the
708    /// block-diagonal-aware fast path is a tracked follow-up.
709    fn jacobian(&self, t: S, z: &[S], jac: &mut [S]) {
710        let n = self.system.n_states();
711        let np = self.system.n_params();
712        let dim = n * (1 + np);
713        let y = &z[..n];
714
715        // Zero the entire augmented matrix first.
716        for slot in jac.iter_mut() {
717            *slot = S::ZERO;
718        }
719
720        let mut jy = self.jy_scratch.borrow_mut();
721        self.system.jacobian_y(t, y, &mut jy);
722
723        // Place J_y on each of the n_s + 1 diagonal blocks.
724        for block in 0..(1 + np) {
725            let row0 = block * n;
726            let col0 = block * n;
727            for i in 0..n {
728                for j in 0..n {
729                    jac[(row0 + i) * dim + (col0 + j)] = jy[i * n + j];
730                }
731            }
732        }
733    }
734
735    // ------------------------------------------------------------------
736    // Mass-matrix + autonomy surface (F-IC-SENS-MASS-SURFACE, 0.1.5,
737    // Foundation Spec §6 #14). The pre-0.1.5 `impl OdeSystem for
738    // AugmentedSystem` inherited the trait defaults for these five
739    // methods — silently dropping the wrapped `ParametricOdeSystem`'s
740    // mass-matrix and autonomy declarations onto the augmented state.
741    // The downstream solvers (Radau5, BDF) read these methods directly
742    // (radau5.rs:234,237; bdf.rs:470,472), so the silent defaults
743    // produced M-blind variational integrations for any host with
744    // non-identity or singular M. The five overrides below delegate /
745    // lift from the host. See `lift_mass_matrix_for_augmented` and
746    // `lift_algebraic_indices_for_augmented` for the block-diagonal
747    // derivation.
748    // ------------------------------------------------------------------
749
750    fn is_autonomous(&self) -> bool {
751        // Variational equations inherit the host's time-dependence
752        // structure — if `f(t, y, p)` is t-independent, so is the
753        // augmented RHS (whose sensitivity columns depend on `J_y`,
754        // `J_p`, themselves evaluated at the host's `(t, y, p)`).
755        self.system.is_autonomous()
756    }
757
758    fn has_mass_matrix(&self) -> bool {
759        // Direct delegation. Downstream solvers (Radau5, BDF) gate the
760        // M-aware Newton system on this flag; a host that declares
761        // `has_mass_matrix = false` must continue to take the identity-M
762        // fast path through the augmented integration.
763        self.system.has_mass_matrix()
764    }
765
766    /// Fill the augmented mass matrix.
767    ///
768    /// Two distinct cases (correctness, not defense): the
769    /// `has_mass_matrix() = false` branch returns the full augmented
770    /// identity rather than the host's identity-lifted result, because
771    /// downstream solvers dispatch on `has_mass_matrix()` separately
772    /// from `mass_matrix()` — a host that explicitly declared
773    /// `has_mass_matrix() = true` with `M = I` is semantically distinct
774    /// from a host that declared `has_mass_matrix() = false`, and the
775    /// branch preserves that distinction. The M-aware branch lifts the
776    /// host's `N × N` `M` block-diagonally over `(N_s + 1)` copies via
777    /// `lift_mass_matrix_for_augmented`.
778    fn mass_matrix(&self, mass: &mut [S]) {
779        let n = self.system.n_states();
780        let np = self.system.n_params();
781        let dim = n * (1 + np);
782        if !self.system.has_mass_matrix() {
783            // Augmented identity for the full `dim × dim` matrix.
784            // Semantically distinct from the M-aware identity branch
785            // (see the rustdoc above).
786            for i in 0..dim {
787                for j in 0..dim {
788                    mass[i * dim + j] = if i == j { S::ONE } else { S::ZERO };
789                }
790            }
791            return;
792        }
793        // Borrow the host's `N × N` M into local scratch and lift onto
794        // the augmented diagonal. Local scratch (not a `RefCell` on the
795        // struct) because this path is cold relative to `rhs`/`jacobian`
796        // — solvers call `mass_matrix` once at the start of integration,
797        // not every step.
798        let mut host_m = vec![S::ZERO; n * n];
799        self.system.mass_matrix(&mut host_m);
800        lift_mass_matrix_for_augmented(&host_m, n, np, mass);
801    }
802
803    fn is_singular_mass(&self) -> bool {
804        // Direct delegation. Singularity lifts through every block (a
805        // zero row of host M is zero in every augmented diagonal block).
806        self.system.is_singular_mass()
807    }
808
809    fn algebraic_indices(&self) -> Vec<usize> {
810        // Cached at construction in `new()`; see the field doc on
811        // `augmented_algebraic_indices` for the host-stability
812        // assumption.
813        self.augmented_algebraic_indices.clone()
814    }
815}
816
817/// Forwarding impl: a reference to a parametric system is itself a parametric
818/// system. Lets [`solve_forward_sensitivity`] accept `&Sys` without taking
819/// ownership of the user's data.
820impl<S: Scalar, T: ParametricOdeSystem<S>> ParametricOdeSystem<S> for &T {
821    fn n_states(&self) -> usize {
822        (*self).n_states()
823    }
824    fn n_params(&self) -> usize {
825        (*self).n_params()
826    }
827    fn params(&self) -> &[S] {
828        (*self).params()
829    }
830    fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]) {
831        (*self).rhs_with_params(t, y, p, dydt)
832    }
833    fn rhs(&self, t: S, y: &[S], dydt: &mut [S]) {
834        (*self).rhs(t, y, dydt)
835    }
836    fn jacobian_y(&self, t: S, y: &[S], jac: &mut [S]) {
837        (*self).jacobian_y(t, y, jac)
838    }
839    fn jacobian_p(&self, t: S, y: &[S], jp: &mut [S]) {
840        (*self).jacobian_p(t, y, jp)
841    }
842    fn initial_sensitivity(&self, y0: &[S], s0: &mut [S]) {
843        (*self).initial_sensitivity(y0, s0)
844    }
845    fn has_analytical_jacobian_y(&self) -> bool {
846        (*self).has_analytical_jacobian_y()
847    }
848    fn has_analytical_jacobian_p(&self) -> bool {
849        (*self).has_analytical_jacobian_p()
850    }
851    // Mass-matrix surface (F-IC-SENS-MASS-SURFACE, 0.1.5). Forwards from
852    // the referenced parametric system; otherwise `&Sys` would silently
853    // inherit the trait defaults and drop the host's M declaration —
854    // exactly the inward-direction boundary gap this PR exists to close
855    // (see `feedback_wrapper_boundary_symmetry`).
856    fn is_autonomous(&self) -> bool {
857        (*self).is_autonomous()
858    }
859    fn has_mass_matrix(&self) -> bool {
860        (*self).has_mass_matrix()
861    }
862    fn mass_matrix(&self, mass: &mut [S]) {
863        (*self).mass_matrix(mass)
864    }
865    fn is_singular_mass(&self) -> bool {
866        (*self).is_singular_mass()
867    }
868    fn algebraic_indices(&self) -> Vec<usize> {
869        (*self).algebraic_indices()
870    }
871}
872
873/// Result of an ODE forward-sensitivity solve.
874///
875/// Layout matches the conventions in the [module-level docs](self):
876/// the state `y` is row-major over time; the sensitivity is row-major over
877/// time, column-major within each per-time block.
878#[derive(Clone, Debug)]
879pub struct SensitivityResult<S: Scalar> {
880    /// Output time points (length `n_times`).
881    pub t: Vec<S>,
882    /// State trajectory: `y[i*N + j] = y_j(t_i)`, length `n_times * N`.
883    pub y: Vec<S>,
884    /// Sensitivities: per-time block of length `N · N_s`, column-major
885    /// within the block, so
886    /// `sensitivity[i * (N*N_s) + k*N + j] = ∂y_j(t_i)/∂p_k`.
887    pub sensitivity: Vec<S>,
888    /// Number of state variables `N`.
889    pub n_states: usize,
890    /// Number of parameters `N_s`.
891    pub n_params: usize,
892    /// Solver statistics from the underlying integration.
893    pub stats: SolverStats,
894    /// Was the integration successful?
895    pub success: bool,
896    /// Error description if the integration failed.
897    pub message: String,
898}
899
900impl<S: Scalar> SensitivityResult<S> {
901    /// Number of output time points.
902    pub fn len(&self) -> usize {
903        self.t.len()
904    }
905
906    /// `true` if there are no output time points.
907    pub fn is_empty(&self) -> bool {
908        self.t.is_empty()
909    }
910
911    /// State vector at output index `i` (length `N`).
912    pub fn y_at(&self, i: usize) -> &[S] {
913        let n = self.n_states;
914        let off = i * n;
915        &self.y[off..off + n]
916    }
917
918    /// Full sensitivity block at output index `i`, column-major
919    /// (length `N · N_s`). `sensitivity_at(i)[k*N + j] = ∂y_j(t_i)/∂p_k`.
920    pub fn sensitivity_at(&self, i: usize) -> &[S] {
921        let block = self.n_states * self.n_params;
922        let off = i * block;
923        &self.sensitivity[off..off + block]
924    }
925
926    /// Contiguous slice `S_{:,k}(t_i) = ∂y/∂p_k` at output index `i`,
927    /// length `N`. Free slice thanks to the column-major flattening.
928    pub fn sensitivity_for_param(&self, i: usize, k: usize) -> &[S] {
929        let n = self.n_states;
930        let block = n * self.n_params;
931        let off = i * block + k * n;
932        &self.sensitivity[off..off + n]
933    }
934
935    /// Single component `∂y_state(t_i)/∂p_param`.
936    pub fn dyi_dpj(&self, i: usize, state: usize, param: usize) -> S {
937        let n = self.n_states;
938        let block = n * self.n_params;
939        self.sensitivity[i * block + param * n + state]
940    }
941
942    /// State at the final output time. Returns an empty slice if the result
943    /// has no time points.
944    pub fn final_state(&self) -> &[S] {
945        if self.is_empty() {
946            &[]
947        } else {
948            self.y_at(self.len() - 1)
949        }
950    }
951
952    /// Sensitivity block at the final output time, column-major (length
953    /// `N · N_s`). Empty slice when the result has no time points.
954    pub fn final_sensitivity(&self) -> &[S] {
955        if self.is_empty() {
956            &[]
957        } else {
958            self.sensitivity_at(self.len() - 1)
959        }
960    }
961
962    /// Normalised (logarithmic) sensitivity at output index `i`:
963    /// `(p_k / y_j) · ∂y_j/∂p_k`, column-major in the same shape as
964    /// [`Self::sensitivity_at`]. Components where `|y_j|` falls below
965    /// `1e-15 · max(1, |y|_∞)` are reported as zero to avoid blow-up.
966    pub fn normalized_sensitivity_at(&self, i: usize, p_nominal: &[S]) -> Vec<S> {
967        assert_eq!(
968            p_nominal.len(),
969            self.n_params,
970            "p_nominal length {} does not match n_params {}",
971            p_nominal.len(),
972            self.n_params,
973        );
974        let n = self.n_states;
975        let np = self.n_params;
976        let y = self.y_at(i);
977        let s = self.sensitivity_at(i);
978
979        let mut y_max = S::ZERO;
980        for &v in y {
981            let av = v.abs();
982            if av > y_max {
983                y_max = av;
984            }
985        }
986        let zero_threshold = S::from_f64(1e-15) * (S::ONE + y_max);
987
988        let mut out = vec![S::ZERO; n * np];
989        for k in 0..np {
990            for j in 0..n {
991                let yj = y[j];
992                if yj.abs() <= zero_threshold {
993                    out[k * n + j] = S::ZERO;
994                } else {
995                    out[k * n + j] = (p_nominal[k] / yj) * s[k * n + j];
996                }
997            }
998        }
999        out
1000    }
1001}
1002
1003/// Closure-shaped adapter for [`solve_forward_sensitivity_with`].
1004///
1005/// Wraps a closure of shape `Fn(t, y, p, dydt)` and a fixed parameter vector
1006/// into a [`ParametricOdeSystem`] using FD-default Jacobians. Public so that
1007/// callers writing a one-shot harness can construct it explicitly when they
1008/// need to share it across multiple solver invocations.
1009pub struct ClosureSystem<S: Scalar, F> {
1010    rhs: F,
1011    params: Vec<S>,
1012    n_states: usize,
1013}
1014
1015impl<S: Scalar, F> ClosureSystem<S, F>
1016where
1017    F: Fn(S, &[S], &[S], &mut [S]),
1018{
1019    /// Build a closure-backed [`ParametricOdeSystem`].
1020    pub fn new(rhs: F, params: Vec<S>, n_states: usize) -> Self {
1021        Self {
1022            rhs,
1023            params,
1024            n_states,
1025        }
1026    }
1027}
1028
1029impl<S: Scalar, F> ParametricOdeSystem<S> for ClosureSystem<S, F>
1030where
1031    F: Fn(S, &[S], &[S], &mut [S]),
1032{
1033    fn n_states(&self) -> usize {
1034        self.n_states
1035    }
1036    fn n_params(&self) -> usize {
1037        self.params.len()
1038    }
1039    fn params(&self) -> &[S] {
1040        &self.params
1041    }
1042    fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]) {
1043        (self.rhs)(t, y, p, dydt)
1044    }
1045}
1046
1047/// Compute the forward sensitivity matrix `S(t) = ∂y(t) / ∂p` of an ODE
1048/// solution using any [`Solver`].
1049///
1050/// Builds the augmented state `z = [y; vec(S)]` (column-major sensitivity
1051/// flattening; see [module-level docs](self)), drives `Sol` over the
1052/// augmented system, and splits the trajectory back into `y` and
1053/// sensitivity blocks at every output time the solver produced.
1054///
1055/// # Errors
1056///
1057/// Returns whatever error `Sol::solve` raises. If the underlying integration
1058/// returns `success = false` (e.g. step-size collapse or max-steps reached),
1059/// the returned [`SensitivityResult`] propagates that with empty trajectory
1060/// vectors and the solver's diagnostic message.
1061///
1062/// # Example
1063///
1064/// ```
1065/// use numra_ode::sensitivity::{solve_forward_sensitivity, ParametricOdeSystem};
1066/// use numra_ode::{DoPri5, SolverOptions};
1067///
1068/// // dy/dt = -k * y, y(0) = 1, k = 0.5.
1069/// // Analytical: y(t) = exp(-k t), ∂y/∂k = -t exp(-k t).
1070/// struct Decay { k: f64 }
1071/// impl ParametricOdeSystem<f64> for Decay {
1072///     fn n_states(&self) -> usize { 1 }
1073///     fn n_params(&self) -> usize { 1 }
1074///     fn params(&self) -> &[f64] { std::slice::from_ref(&self.k) }
1075///     fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
1076///         dy[0] = -p[0] * y[0];
1077///     }
1078///     fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) { jy[0] = -self.k; }
1079///     fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) { jp[0] = -y[0]; }
1080///     fn has_analytical_jacobian_y(&self) -> bool { true }
1081///     fn has_analytical_jacobian_p(&self) -> bool { true }
1082/// }
1083///
1084/// let r = solve_forward_sensitivity::<DoPri5, _, _>(
1085///     &Decay { k: 0.5 },
1086///     0.0, 2.0,
1087///     &[1.0],
1088///     &SolverOptions::default().rtol(1e-8).atol(1e-10),
1089/// ).unwrap();
1090///
1091/// let last = r.len() - 1;
1092/// let dy_dk = r.dyi_dpj(last, 0, 0);
1093/// let analytical = -r.t[last] * (-0.5 * r.t[last]).exp();
1094/// assert!((dy_dk - analytical).abs() < 1e-5);
1095/// ```
1096pub fn solve_forward_sensitivity<Sol, S, Sys>(
1097    system: &Sys,
1098    t0: S,
1099    tf: S,
1100    y0: &[S],
1101    options: &SolverOptions<S>,
1102) -> Result<SensitivityResult<S>, SolverError>
1103where
1104    Sol: Solver<S>,
1105    S: Scalar,
1106    Sys: ParametricOdeSystem<S>,
1107{
1108    let n = system.n_states();
1109    let np = system.n_params();
1110    assert_eq!(
1111        y0.len(),
1112        n,
1113        "solve_forward_sensitivity: y0.len() = {} but n_states = {}",
1114        y0.len(),
1115        n,
1116    );
1117    assert_eq!(
1118        system.params().len(),
1119        np,
1120        "solve_forward_sensitivity: params().len() = {} but n_params = {}",
1121        system.params().len(),
1122        np,
1123    );
1124
1125    let aug = AugmentedSystem::new(system);
1126    let z0 = aug.initial_augmented(y0);
1127    let aug_dim = aug.augmented_dim();
1128
1129    let aug_result = Sol::solve(&aug, t0, tf, &z0, options)?;
1130
1131    if !aug_result.success {
1132        return Ok(SensitivityResult {
1133            t: Vec::new(),
1134            y: Vec::new(),
1135            sensitivity: Vec::new(),
1136            n_states: n,
1137            n_params: np,
1138            stats: aug_result.stats,
1139            success: false,
1140            message: aug_result.message,
1141        });
1142    }
1143
1144    let n_times = aug_result.len();
1145    let mut t = Vec::with_capacity(n_times);
1146    let mut y = Vec::with_capacity(n_times * n);
1147    let mut sensitivity = Vec::with_capacity(n_times * n * np);
1148
1149    for i in 0..n_times {
1150        t.push(aug_result.t[i]);
1151        let block_start = i * aug_dim;
1152        let state_block = &aug_result.y[block_start..block_start + n];
1153        let sens_block = &aug_result.y[block_start + n..block_start + aug_dim];
1154        y.extend_from_slice(state_block);
1155        sensitivity.extend_from_slice(sens_block);
1156    }
1157
1158    Ok(SensitivityResult {
1159        t,
1160        y,
1161        sensitivity,
1162        n_states: n,
1163        n_params: np,
1164        stats: aug_result.stats,
1165        success: true,
1166        message: String::new(),
1167    })
1168}
1169
1170/// Closure-shaped convenience wrapper around [`solve_forward_sensitivity`].
1171///
1172/// Wraps `rhs(t, y, p, dydt)` plus a parameter vector into an internal
1173/// [`ClosureSystem`] (FD-default Jacobians) and forwards. Suitable for
1174/// one-shot analyses, tests, and REPL-style scripts. For production usage
1175/// — particularly stiff problems where analytical Jacobians materially
1176/// change runtime — implement [`ParametricOdeSystem`] directly and call
1177/// [`solve_forward_sensitivity`].
1178///
1179/// # Example
1180///
1181/// ```
1182/// use numra_ode::{DoPri5, SolverOptions};
1183/// use numra_ode::sensitivity::solve_forward_sensitivity_with;
1184///
1185/// // dy/dt = -k * y, y(0) = 1, k = 0.5.
1186/// let r = solve_forward_sensitivity_with::<DoPri5, f64, _>(
1187///     |_t, y, p, dy| { dy[0] = -p[0] * y[0]; },
1188///     &[1.0],
1189///     &[0.5],
1190///     0.0, 2.0,
1191///     &SolverOptions::default().rtol(1e-8).atol(1e-10),
1192/// ).unwrap();
1193///
1194/// let last = r.len() - 1;
1195/// let dy_dk = r.dyi_dpj(last, 0, 0);
1196/// let analytical: f64 = -r.t[last] * (-0.5 * r.t[last]).exp();
1197/// assert!((dy_dk - analytical).abs() < 1e-4);
1198/// ```
1199pub fn solve_forward_sensitivity_with<Sol, S, F>(
1200    rhs: F,
1201    y0: &[S],
1202    params: &[S],
1203    t0: S,
1204    tf: S,
1205    options: &SolverOptions<S>,
1206) -> Result<SensitivityResult<S>, SolverError>
1207where
1208    Sol: Solver<S>,
1209    S: Scalar,
1210    F: Fn(S, &[S], &[S], &mut [S]),
1211{
1212    let system = ClosureSystem {
1213        rhs,
1214        params: params.to_vec(),
1215        n_states: y0.len(),
1216    };
1217    solve_forward_sensitivity::<Sol, S, _>(&system, t0, tf, y0, options)
1218}
1219
1220// ============================================================================
1221// Initial-condition sensitivity (state-transition matrix)
1222//
1223// Specialises the variational machinery above to the IC case: integrates
1224// `Ṡ = J_y · S` with `S(t₀) = I_N`, producing `Φ(t) = ∂y(t)/∂y₀` —
1225// the state-transition matrix of the linearised flow. Reuses
1226// `AugmentedSystem` and `solve_forward_sensitivity` unchanged; no fork.
1227//
1228// The user-facing API takes an `OdeSystem<S>` (or RHS closure) and returns
1229// a `StateTransitionResult<S>` whose accessors speak Φ(t)-vocabulary only.
1230// Internally a private `IcAsParametric` wrapper reinterprets the plain
1231// system as a `ParametricOdeSystem` with `n_params := n_states`,
1232// `J_p ≡ 0`, and `S(t₀) = I`. The required dummy `params()` slice is the
1233// only seam, encapsulated behind the wrapper and pinned by the test
1234// `ic_dummy_params_are_unread` below.
1235//
1236// Invariant **IC-NO-PARAM-READ.** `IcAsParametric::rhs_with_params(t,y,p,dy)`
1237// never reads `p`. The wrapped `OdeSystem::rhs(t,y,dy)` has no `p` argument
1238// and cannot. Therefore `∂f/∂p ≡ 0` is mathematically exact (not "small in
1239// some norm"), and `IcAsParametric::jacobian_p` writing zeros with
1240// `has_analytical_jacobian_p() = true` is a correct analytical declaration
1241// — `AugmentedSystem`'s debug consistency check at `sensitivity.rs:373` is
1242// truthfully suppressed by this flag, not silenced.
1243// ============================================================================
1244
1245/// Private wrapper: makes an `OdeSystem<S>` look like a `ParametricOdeSystem<S>`
1246/// whose "parameters" are seeded so the forward-sensitivity flow computes
1247/// `Φ(t) = ∂y(t)/∂y₀` instead of `∂y(t)/∂p`.
1248struct IcAsParametric<'a, S: Scalar, Sys: OdeSystem<S> + ?Sized> {
1249    inner: &'a Sys,
1250    /// Dummy parameter slice of length `n_states`. Length-only requirement
1251    /// from `ParametricOdeSystem` (asserted at solve entry, line 852–858).
1252    /// Provably unread by `rhs_with_params`; pinned by
1253    /// `ic_dummy_params_are_unread`.
1254    dummy_p: Vec<S>,
1255}
1256
1257impl<'a, S: Scalar, Sys: OdeSystem<S> + ?Sized> ParametricOdeSystem<S>
1258    for IcAsParametric<'a, S, Sys>
1259{
1260    fn n_states(&self) -> usize {
1261        self.inner.dim()
1262    }
1263    fn n_params(&self) -> usize {
1264        self.inner.dim()
1265    }
1266    fn params(&self) -> &[S] {
1267        &self.dummy_p
1268    }
1269
1270    fn rhs_with_params(&self, t: S, y: &[S], _p: &[S], dy: &mut [S]) {
1271        // IC-NO-PARAM-READ: `_p` is provably unread; `inner.rhs` has no
1272        // `p` argument and cannot consume it.
1273        self.inner.rhs(t, y, dy);
1274    }
1275
1276    fn jacobian_y(&self, t: S, y: &[S], jy: &mut [S]) {
1277        // Forwards to the inner system's Jacobian — analytical override if
1278        // the user supplied one (e.g. `AutodiffJacobianSystem` from the
1279        // `autodiff` feature), or `OdeSystem`'s FD default otherwise.
1280        self.inner.jacobian(t, y, jy);
1281    }
1282
1283    fn jacobian_p(&self, _t: S, _y: &[S], jp: &mut [S]) {
1284        // J_p ≡ 0 analytically (IC-NO-PARAM-READ ⇒ ∂f/∂p ≡ 0 exactly).
1285        for slot in jp.iter_mut() {
1286            *slot = S::ZERO;
1287        }
1288    }
1289
1290    fn initial_sensitivity(&self, _y0: &[S], s0: &mut [S]) {
1291        // Column-major identity: s0[k*N + i] = δ_{ik}.
1292        let n = self.inner.dim();
1293        for slot in s0.iter_mut() {
1294            *slot = S::ZERO;
1295        }
1296        for k in 0..n {
1297            s0[k * n + k] = S::ONE;
1298        }
1299    }
1300
1301    fn has_analytical_jacobian_y(&self) -> bool {
1302        // Claim analytical so AugmentedSystem honours the forwarded
1303        // `inner.jacobian` (which may itself be analytical or FD-default).
1304        // The debug consistency check is then irrelevant for this flag.
1305        true
1306    }
1307
1308    fn has_analytical_jacobian_p(&self) -> bool {
1309        // Truthful: J_p ≡ 0 is analytically exact per IC-NO-PARAM-READ.
1310        true
1311    }
1312
1313    // ------------------------------------------------------------------
1314    // Mass-matrix + autonomy forwarding (F-IC-SENS-MASS-SURFACE, 0.1.5,
1315    // Foundation Spec §6 #14). Closes the inward-direction boundary gap
1316    // F-IC-SENS shipped with: the 0.1.4 wrapper inherited the
1317    // `ParametricOdeSystem` defaults for these five methods (returning
1318    // `has_mass_matrix = false`, identity M, etc.), silently dropping
1319    // the wrapped `OdeSystem`'s declarations on the way into the
1320    // augmented integration. See `feedback_wrapper_boundary_symmetry`
1321    // for the both-direction boundary discipline this restores.
1322    // ------------------------------------------------------------------
1323
1324    fn is_autonomous(&self) -> bool {
1325        self.inner.is_autonomous()
1326    }
1327
1328    fn has_mass_matrix(&self) -> bool {
1329        self.inner.has_mass_matrix()
1330    }
1331
1332    fn mass_matrix(&self, mass: &mut [S]) {
1333        self.inner.mass_matrix(mass);
1334    }
1335
1336    fn is_singular_mass(&self) -> bool {
1337        self.inner.is_singular_mass()
1338    }
1339
1340    fn algebraic_indices(&self) -> Vec<usize> {
1341        self.inner.algebraic_indices()
1342    }
1343}
1344
1345/// Result of an initial-condition sensitivity solve.
1346///
1347/// Carries the state trajectory `y(t)` and the state-transition matrix
1348/// `Φ(t) = ∂y(t) / ∂y₀ ∈ ℝ^{N×N}` of the linearised flow, evaluated at every
1349/// output time the underlying solver produced.
1350///
1351/// # Φ(t) layout
1352///
1353/// Column-major within each per-time block: `phi(i)[k*N + j] = Φ(t_i)_{j,k}
1354/// = ∂y_j(t_i) / ∂y_{0,k}`. The column-major flattening makes the j-th
1355/// column (the trajectory's response at time `t_i` to a unit perturbation
1356/// of `y_{0,k}`) a free contiguous slice — see [`Self::phi_column`].
1357///
1358/// # Vocabulary
1359///
1360/// All accessors speak state-transition-matrix vocabulary
1361/// (`Φ`, `∂y_i/∂y_{0,j}`, monodromy). The type holds the underlying
1362/// [`SensitivityResult`] in a private field; it does **not** `Deref` to it
1363/// and does **not** re-export the parameter-shaped accessors.
1364#[derive(Clone, Debug)]
1365pub struct StateTransitionResult<S: Scalar> {
1366    inner: SensitivityResult<S>,
1367}
1368
1369impl<S: Scalar> StateTransitionResult<S> {
1370    /// Output time points (length `n_times`). `t()[i]` is the i-th output
1371    /// time; `phi(i)` is `Φ(t()[i])`.
1372    pub fn t(&self) -> &[S] {
1373        &self.inner.t
1374    }
1375
1376    /// Number of output time points.
1377    pub fn len(&self) -> usize {
1378        self.inner.len()
1379    }
1380
1381    /// `true` iff the result has no output time points.
1382    pub fn is_empty(&self) -> bool {
1383        self.inner.is_empty()
1384    }
1385
1386    /// Underlying state dimension `N`.
1387    pub fn dim(&self) -> usize {
1388        self.inner.n_states
1389    }
1390
1391    /// State `y(t_i)` at output index `i` (length `N`).
1392    pub fn y(&self, i: usize) -> &[S] {
1393        self.inner.y_at(i)
1394    }
1395
1396    /// State at the final output time, length `N`. Empty slice when the
1397    /// result has no time points.
1398    pub fn final_y(&self) -> &[S] {
1399        self.inner.final_state()
1400    }
1401
1402    /// Full state-transition matrix `Φ(t_i) ∈ ℝ^{N×N}` at output index `i`,
1403    /// column-major. `phi(i)[k*N + j] = Φ(t_i)_{j,k} = ∂y_j(t_i)/∂y_{0,k}`.
1404    pub fn phi(&self, i: usize) -> &[S] {
1405        self.inner.sensitivity_at(i)
1406    }
1407
1408    /// Single entry `Φ(t_i)_{row, col} = ∂y_row(t_i) / ∂y_{0, col}`.
1409    pub fn phi_ij(&self, i: usize, row: usize, col: usize) -> S {
1410        self.inner.dyi_dpj(i, row, col)
1411    }
1412
1413    /// Contiguous slice `Φ(t_i)_{:, col}` of length `N` — the trajectory's
1414    /// response at `t_i` to a unit perturbation of `y_{0, col}`. Free slice
1415    /// thanks to the column-major flattening.
1416    pub fn phi_column(&self, i: usize, col: usize) -> &[S] {
1417        self.inner.sensitivity_for_param(i, col)
1418    }
1419
1420    /// `Φ(t_f)`, the state-transition matrix at the final output time.
1421    /// When `(t_f − t_0)` equals one period of a periodic orbit this is
1422    /// the *monodromy matrix*, whose eigenvalues are the orbit's
1423    /// characteristic (Floquet) multipliers. The type does not enforce
1424    /// periodicity; it only names the natural use. Empty slice when the
1425    /// result has no time points.
1426    pub fn final_phi(&self) -> &[S] {
1427        self.inner.final_sensitivity()
1428    }
1429
1430    /// `true` iff the underlying integration succeeded.
1431    pub fn success(&self) -> bool {
1432        self.inner.success
1433    }
1434
1435    /// Solver diagnostic message (empty on success).
1436    pub fn message(&self) -> &str {
1437        &self.inner.message
1438    }
1439
1440    /// Solver statistics from the underlying integration.
1441    pub fn stats(&self) -> &SolverStats {
1442        &self.inner.stats
1443    }
1444}
1445
1446/// Compute the state-transition matrix `Φ(t) = ∂y(t) / ∂y₀` of the
1447/// linearised flow of `dy/dt = f(t, y)`, using any [`Solver`].
1448///
1449/// Integrates the variational equations `Ṡ = J_y · S` with seed
1450/// `S(t₀) = I_N` over the same [`AugmentedSystem`] used for parameter
1451/// sensitivity, then extracts the state-transition trajectory.
1452///
1453/// The user supplies only the ODE system; the IC seeding (identity initial
1454/// sensitivity, zero parameter forcing) is performed internally. `J_y` is
1455/// obtained from [`OdeSystem::jacobian`] — the user's analytical override
1456/// if any, the trait's forward-FD default otherwise. For an
1457/// exact-to-round-off `J_y` from autodiff, wrap the system in
1458/// `AutodiffJacobianSystem` (behind the `autodiff` feature).
1459///
1460/// # Example: scalar exponential decay
1461///
1462/// `dy/dt = -k y` has closed-form `Φ(t) = exp(-k t)`.
1463///
1464/// ```
1465/// use numra_ode::sensitivity::solve_initial_condition_sensitivity_with;
1466/// use numra_ode::{DoPri5, SolverOptions};
1467///
1468/// let k = 0.5_f64;
1469/// let result = solve_initial_condition_sensitivity_with::<DoPri5, f64, _>(
1470///     move |_t, y, dy| { dy[0] = -k * y[0]; },
1471///     &[1.0],
1472///     0.0, 2.0,
1473///     &SolverOptions::default().rtol(1e-9).atol(1e-12),
1474/// ).unwrap();
1475///
1476/// let last = result.len() - 1;
1477/// let phi = result.phi_ij(last, 0, 0);
1478/// let analytical = (-k * result.t()[last]).exp();
1479/// assert!((phi - analytical).abs() < 1e-7);
1480/// ```
1481pub fn solve_initial_condition_sensitivity<Sol, S, Sys>(
1482    system: &Sys,
1483    t0: S,
1484    tf: S,
1485    y0: &[S],
1486    options: &SolverOptions<S>,
1487) -> Result<StateTransitionResult<S>, SolverError>
1488where
1489    Sol: Solver<S>,
1490    S: Scalar,
1491    Sys: OdeSystem<S> + ?Sized,
1492{
1493    let n = system.dim();
1494    assert_eq!(
1495        y0.len(),
1496        n,
1497        "solve_initial_condition_sensitivity: y0.len() = {} but dim() = {}",
1498        y0.len(),
1499        n,
1500    );
1501
1502    let wrapper = IcAsParametric {
1503        inner: system,
1504        dummy_p: vec![S::ZERO; n],
1505    };
1506    let inner = solve_forward_sensitivity::<Sol, S, _>(&wrapper, t0, tf, y0, options)?;
1507    Ok(StateTransitionResult { inner })
1508}
1509
1510/// Closure-shaped convenience wrapper around
1511/// [`solve_initial_condition_sensitivity`].
1512///
1513/// Wraps `rhs(t, y, dydt)` into an internal [`OdeProblem`] (forward-FD
1514/// Jacobian via the [`OdeSystem`] default) and forwards. Suitable for
1515/// one-shot analyses and REPL-style scripts. For an analytical `J_y` or an
1516/// autodiff-derived `J_y`, implement [`OdeSystem`] directly (with a
1517/// `jacobian` override or by wrapping in `AutodiffJacobianSystem`)
1518/// and call [`solve_initial_condition_sensitivity`].
1519pub fn solve_initial_condition_sensitivity_with<Sol, S, F>(
1520    rhs: F,
1521    y0: &[S],
1522    t0: S,
1523    tf: S,
1524    options: &SolverOptions<S>,
1525) -> Result<StateTransitionResult<S>, SolverError>
1526where
1527    Sol: Solver<S>,
1528    S: Scalar,
1529    F: Fn(S, &[S], &mut [S]),
1530{
1531    let problem = OdeProblem::new(rhs, t0, tf, y0.to_vec());
1532    solve_initial_condition_sensitivity::<Sol, S, _>(&problem, t0, tf, y0, options)
1533}
1534
1535#[cfg(test)]
1536mod tests {
1537    use super::*;
1538
1539    /// Minimal scalar decay system used to drive the trait through
1540    /// `OdeSystem::rhs` and `OdeSystem::jacobian` without involving a solver.
1541    struct ExpDecay {
1542        k: f64,
1543    }
1544
1545    impl ParametricOdeSystem<f64> for ExpDecay {
1546        fn n_states(&self) -> usize {
1547            1
1548        }
1549        fn n_params(&self) -> usize {
1550            1
1551        }
1552        fn params(&self) -> &[f64] {
1553            std::slice::from_ref(&self.k)
1554        }
1555        fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
1556            dy[0] = -p[0] * y[0];
1557        }
1558        fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) {
1559            jy[0] = -self.k;
1560        }
1561        fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) {
1562            jp[0] = -y[0];
1563        }
1564        fn has_analytical_jacobian_y(&self) -> bool {
1565            true
1566        }
1567        fn has_analytical_jacobian_p(&self) -> bool {
1568            true
1569        }
1570    }
1571
1572    /// 2-state, 2-parameter linear system (mirrors the module-doc example).
1573    /// Lets us exercise the asymmetric J_y / J_p layouts in code.
1574    struct Lin2 {
1575        p: [f64; 2],
1576    }
1577
1578    impl ParametricOdeSystem<f64> for Lin2 {
1579        fn n_states(&self) -> usize {
1580            2
1581        }
1582        fn n_params(&self) -> usize {
1583            2
1584        }
1585        fn params(&self) -> &[f64] {
1586            &self.p
1587        }
1588        fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
1589            dy[0] = -p[0] * y[0] + p[1] * y[1];
1590            dy[1] = -p[1] * y[1];
1591        }
1592        fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) {
1593            jy[0 * 2 + 0] = -self.p[0];
1594            jy[0 * 2 + 1] = self.p[1];
1595            jy[1 * 2 + 0] = 0.0;
1596            jy[1 * 2 + 1] = -self.p[1];
1597        }
1598        fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) {
1599            jp[0 * 2 + 0] = -y[0];
1600            jp[0 * 2 + 1] = 0.0;
1601            jp[1 * 2 + 0] = y[1];
1602            jp[1 * 2 + 1] = -y[1];
1603        }
1604        fn has_analytical_jacobian_y(&self) -> bool {
1605            true
1606        }
1607        fn has_analytical_jacobian_p(&self) -> bool {
1608            true
1609        }
1610    }
1611
1612    #[test]
1613    fn augmented_dim_and_initial_state() {
1614        let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
1615        assert_eq!(aug.augmented_dim(), 2 * (1 + 2));
1616        let z0 = aug.initial_augmented(&[1.0, 2.0]);
1617        assert_eq!(z0.len(), 6);
1618        assert_eq!(&z0[..2], &[1.0, 2.0]);
1619        for s in &z0[2..] {
1620            assert_eq!(*s, 0.0);
1621        }
1622    }
1623
1624    #[test]
1625    fn augmented_rhs_matches_hand_derivation() {
1626        // System: dy/dt = -k y, S_{0,0} = ∂y/∂k.
1627        // dS/dt = J_y * S + J_p = -k * S - y.
1628        let aug = AugmentedSystem::new(ExpDecay { k: 0.5 });
1629        let z = vec![2.0, 1.0]; // y = 2.0, S = 1.0
1630        let mut dz = vec![0.0; 2];
1631        aug.rhs(0.0, &z, &mut dz);
1632
1633        // dy/dt = -0.5 * 2.0 = -1.0
1634        assert!((dz[0] + 1.0).abs() < 1e-12);
1635        // dS/dt = -0.5 * 1.0 + (-2.0) = -2.5
1636        assert!((dz[1] + 2.5).abs() < 1e-12);
1637    }
1638
1639    #[test]
1640    fn augmented_rhs_matches_hand_derivation_two_param() {
1641        // Lin2 system at y = (1, 2), p = (1, 0.5), S = identity (initial sensitivity ID).
1642        let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
1643
1644        // z = [y_0, y_1, S_{0,0}, S_{1,0}, S_{0,1}, S_{1,1}]
1645        //     [ 1,    2,    1,       0,       0,       1   ]
1646        let z = vec![1.0, 2.0, 1.0, 0.0, 0.0, 1.0];
1647        let mut dz = vec![0.0; 6];
1648        aug.rhs(0.0, &z, &mut dz);
1649
1650        // dy/dt
1651        // dy_0/dt = -1*1 + 0.5*2 = 0.0
1652        // dy_1/dt =        -0.5*2 = -1.0
1653        assert!((dz[0] - 0.0).abs() < 1e-12);
1654        assert!((dz[1] + 1.0).abs() < 1e-12);
1655
1656        // dS_{:,0}/dt = J_y · [1; 0] + J_p_{:,0}
1657        // J_y · [1; 0] = [-1; 0];   J_p_{:,0} = [-1; 0]; total = [-2; 0]
1658        assert!((dz[2] + 2.0).abs() < 1e-12);
1659        assert!((dz[3] - 0.0).abs() < 1e-12);
1660
1661        // dS_{:,1}/dt = J_y · [0; 1] + J_p_{:,1}
1662        // J_y · [0; 1] = [0.5; -0.5]; J_p_{:,1} = [2; -2]; total = [2.5; -2.5]
1663        assert!((dz[4] - 2.5).abs() < 1e-12);
1664        assert!((dz[5] + 2.5).abs() < 1e-12);
1665    }
1666
1667    #[test]
1668    fn augmented_jacobian_is_block_diagonal_with_jy() {
1669        let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
1670        let z = vec![1.0, 2.0, 1.0, 0.0, 0.0, 1.0];
1671        let dim = aug.augmented_dim();
1672        let mut jac = vec![0.0; dim * dim];
1673        aug.jacobian(0.0, &z, &mut jac);
1674
1675        // Expected J_y at p=(1, 0.5):
1676        //   [-1.0  0.5]
1677        //   [ 0.0 -0.5]
1678        let jy = [-1.0_f64, 0.5, 0.0, -0.5];
1679
1680        for block in 0..3 {
1681            let row0 = block * 2;
1682            let col0 = block * 2;
1683            for i in 0..2 {
1684                for j in 0..2 {
1685                    let got = jac[(row0 + i) * dim + (col0 + j)];
1686                    let want = jy[i * 2 + j];
1687                    assert!(
1688                        (got - want).abs() < 1e-12,
1689                        "diag block {block}, ({i},{j}): got {got}, want {want}"
1690                    );
1691                }
1692            }
1693        }
1694
1695        // Off-diagonal blocks must be zero.
1696        for br in 0..3 {
1697            for bc in 0..3 {
1698                if br == bc {
1699                    continue;
1700                }
1701                for i in 0..2 {
1702                    for j in 0..2 {
1703                        let r = br * 2 + i;
1704                        let c = bc * 2 + j;
1705                        assert_eq!(jac[r * dim + c], 0.0, "off-diag ({br},{bc}) ({i},{j})");
1706                    }
1707                }
1708            }
1709        }
1710    }
1711
1712    #[test]
1713    fn fd_default_jacobian_y_matches_analytical() {
1714        // Re-implement Lin2 without analytical J_y so the FD default is exercised.
1715        struct Lin2Fd {
1716            p: [f64; 2],
1717        }
1718        impl ParametricOdeSystem<f64> for Lin2Fd {
1719            fn n_states(&self) -> usize {
1720                2
1721            }
1722            fn n_params(&self) -> usize {
1723                2
1724            }
1725            fn params(&self) -> &[f64] {
1726                &self.p
1727            }
1728            fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
1729                dy[0] = -p[0] * y[0] + p[1] * y[1];
1730                dy[1] = -p[1] * y[1];
1731            }
1732            // No jacobian_y / jacobian_p override — FD defaults are used.
1733        }
1734
1735        let analytical = Lin2 { p: [1.0, 0.5] };
1736        let fd = Lin2Fd { p: [1.0, 0.5] };
1737
1738        let y = [1.0, 2.0];
1739        let mut a_jy = [0.0; 4];
1740        let mut f_jy = [0.0; 4];
1741        analytical.jacobian_y(0.0, &y, &mut a_jy);
1742        fd.jacobian_y(0.0, &y, &mut f_jy);
1743        for k in 0..4 {
1744            assert!(
1745                (a_jy[k] - f_jy[k]).abs() < 1e-7,
1746                "j_y[{k}] analytical={} fd={}",
1747                a_jy[k],
1748                f_jy[k]
1749            );
1750        }
1751
1752        let mut a_jp = [0.0; 4];
1753        let mut f_jp = [0.0; 4];
1754        analytical.jacobian_p(0.0, &y, &mut a_jp);
1755        fd.jacobian_p(0.0, &y, &mut f_jp);
1756        for k in 0..4 {
1757            assert!(
1758                (a_jp[k] - f_jp[k]).abs() < 1e-7,
1759                "j_p[{k}] analytical={} fd={}",
1760                a_jp[k],
1761                f_jp[k]
1762            );
1763        }
1764    }
1765
1766    #[test]
1767    fn sensitivity_result_accessors() {
1768        // Hand-construct a result for a 2-state, 3-param trajectory at 2 time points.
1769        // sensitivity layout per time block (column-major, length N*N_s = 6):
1770        //   [ S_{0,0}  S_{1,0}  | S_{0,1}  S_{1,1}  | S_{0,2}  S_{1,2} ]
1771        let res = SensitivityResult {
1772            t: vec![0.0, 1.0],
1773            y: vec![1.0, 2.0, 0.5, 1.5],
1774            sensitivity: vec![
1775                // t = 0
1776                0.1, 0.2, 0.3, 0.4, 0.5, 0.6, // t = 1
1777                1.1, 1.2, 1.3, 1.4, 1.5, 1.6,
1778            ],
1779            n_states: 2,
1780            n_params: 3,
1781            stats: SolverStats::new(),
1782            success: true,
1783            message: String::new(),
1784        };
1785
1786        assert_eq!(res.len(), 2);
1787        assert_eq!(res.y_at(1), &[0.5, 1.5]);
1788        assert_eq!(res.sensitivity_at(0).len(), 6);
1789        assert_eq!(res.sensitivity_for_param(1, 2), &[1.5, 1.6]);
1790        assert!((res.dyi_dpj(0, 1, 0) - 0.2).abs() < 1e-12);
1791        assert!((res.dyi_dpj(1, 0, 2) - 1.5).abs() < 1e-12);
1792        assert_eq!(res.final_state(), &[0.5, 1.5]);
1793        assert_eq!(res.final_sensitivity().len(), 6);
1794    }
1795
1796    #[test]
1797    fn normalized_sensitivity_handles_zero_state_safely() {
1798        let res = SensitivityResult {
1799            t: vec![0.0],
1800            y: vec![1.0, 0.0],
1801            sensitivity: vec![0.5, 0.7],
1802            n_states: 2,
1803            n_params: 1,
1804            stats: SolverStats::new(),
1805            success: true,
1806            message: String::new(),
1807        };
1808        let norm = res.normalized_sensitivity_at(0, &[2.0]);
1809        // (p / y_0) · S_{0,0} = (2 / 1) · 0.5 = 1.0
1810        assert!((norm[0] - 1.0).abs() < 1e-12);
1811        // y_1 = 0 → reported as zero rather than ∞.
1812        assert_eq!(norm[1], 0.0);
1813    }
1814
1815    // ---------------------------------------------------------------
1816    // Solve entry-point sanity tests.
1817    //
1818    // The full solver-coverage matrix lives in
1819    // numra-ode/tests/sensitivity_regression.rs (commit 3).
1820    // ---------------------------------------------------------------
1821
1822    use crate::DoPri5;
1823    use crate::SolverOptions;
1824
1825    #[test]
1826    fn solve_forward_sensitivity_with_dopri5_matches_analytical_decay() {
1827        // dy/dt = -k y, y(0) = 1, k = 0.5.
1828        // Analytical: ∂y/∂k(t) = -t · exp(-k t).
1829        let result = solve_forward_sensitivity_with::<DoPri5, f64, _>(
1830            |_t, y, p, dy| {
1831                dy[0] = -p[0] * y[0];
1832            },
1833            &[1.0],
1834            &[0.5],
1835            0.0,
1836            2.0,
1837            &SolverOptions::default().rtol(1e-8).atol(1e-10),
1838        )
1839        .expect("solve failed");
1840
1841        assert!(result.success);
1842        assert_eq!(result.n_states, 1);
1843        assert_eq!(result.n_params, 1);
1844        let last = result.len() - 1;
1845        let dy_dk = result.dyi_dpj(last, 0, 0);
1846        let t_last = result.t[last];
1847        let analytical = -t_last * (-0.5 * t_last).exp();
1848        assert!(
1849            (dy_dk - analytical).abs() < 1e-5,
1850            "computed {dy_dk}, analytical {analytical}, |err| = {}",
1851            (dy_dk - analytical).abs()
1852        );
1853    }
1854
1855    #[test]
1856    fn solve_trait_form_matches_closure_form() {
1857        let trait_result = solve_forward_sensitivity::<DoPri5, _, _>(
1858            &ExpDecay { k: 0.5 },
1859            0.0,
1860            2.0,
1861            &[1.0],
1862            &SolverOptions::default().rtol(1e-9).atol(1e-12),
1863        )
1864        .expect("trait solve failed");
1865
1866        let closure_result = solve_forward_sensitivity_with::<DoPri5, f64, _>(
1867            |_t, y, p, dy| {
1868                dy[0] = -p[0] * y[0];
1869            },
1870            &[1.0],
1871            &[0.5],
1872            0.0,
1873            2.0,
1874            &SolverOptions::default().rtol(1e-9).atol(1e-12),
1875        )
1876        .expect("closure solve failed");
1877
1878        // Trait form has analytical Jacobians; closure form falls back to FD.
1879        // Both must converge to the same answer at the final time, modulo FD
1880        // noise.
1881        let n_t = trait_result.len();
1882        let n_c = closure_result.len();
1883        let trait_final = trait_result.dyi_dpj(n_t - 1, 0, 0);
1884        let closure_final = closure_result.dyi_dpj(n_c - 1, 0, 0);
1885        assert!(
1886            (trait_final - closure_final).abs() < 1e-5,
1887            "trait {trait_final} vs closure {closure_final}, |err| = {}",
1888            (trait_final - closure_final).abs()
1889        );
1890    }
1891
1892    #[test]
1893    fn solve_propagates_solver_failure_as_unsuccessful_result() {
1894        // Forcing max_steps = 1 against a multi-step problem makes the
1895        // underlying solver return success = false (or an error). We assert
1896        // that solve_forward_sensitivity surfaces this without panicking.
1897        let opts = SolverOptions::default().rtol(1e-9).atol(1e-12).max_steps(1);
1898
1899        let outcome = solve_forward_sensitivity_with::<DoPri5, f64, _>(
1900            |_t, y, p, dy| {
1901                dy[0] = -p[0] * y[0];
1902            },
1903            &[1.0],
1904            &[0.5],
1905            0.0,
1906            10.0,
1907            &opts,
1908        );
1909
1910        match outcome {
1911            Ok(r) => {
1912                assert!(!r.success, "expected unsuccessful result, got success");
1913                assert!(r.t.is_empty());
1914                assert!(r.y.is_empty());
1915                assert!(r.sensitivity.is_empty());
1916            }
1917            Err(_) => {
1918                // Some solvers return Err directly for max_steps; that's also acceptable.
1919            }
1920        }
1921    }
1922
1923    // ====================================================================
1924    // F-IC-SENS-MASS-SURFACE Test 6 — structural: IcAsParametric forwards
1925    // the wrapped OdeSystem's mass-matrix surface and is_autonomous.
1926    //
1927    // Closes the inward-direction structural check of the boundary-
1928    // symmetry pair. Lives here (not in tests/) because IcAsParametric
1929    // is private. Paired with Test 5 in
1930    // `numra-ode/tests/ic_sensitivity_mass_matrix.rs::augmented_system_
1931    // lifts_mass_surface_block_diagonally`, which covers the
1932    // AugmentedSystem (outward-from-wrapper) lift on the same data.
1933    // ====================================================================
1934
1935    struct MockMassOde;
1936
1937    impl OdeSystem<f64> for MockMassOde {
1938        fn dim(&self) -> usize {
1939            3
1940        }
1941        fn rhs(&self, _t: f64, _y: &[f64], dy: &mut [f64]) {
1942            for slot in dy.iter_mut() {
1943                *slot = 0.0;
1944            }
1945        }
1946        fn is_autonomous(&self) -> bool {
1947            true
1948        }
1949        fn has_mass_matrix(&self) -> bool {
1950            true
1951        }
1952        fn mass_matrix(&self, m: &mut [f64]) {
1953            // diag(2, 3, 0) row-major.
1954            for slot in m.iter_mut() {
1955                *slot = 0.0;
1956            }
1957            m[0 * 3 + 0] = 2.0;
1958            m[1 * 3 + 1] = 3.0;
1959            m[2 * 3 + 2] = 0.0;
1960        }
1961        fn is_singular_mass(&self) -> bool {
1962            true
1963        }
1964        fn algebraic_indices(&self) -> Vec<usize> {
1965            vec![2]
1966        }
1967    }
1968
1969    #[test]
1970    fn ic_as_parametric_forwards_mass_surface() {
1971        let host = MockMassOde;
1972        let n = host.dim();
1973        let wrapper = IcAsParametric {
1974            inner: &host,
1975            dummy_p: vec![0.0; n],
1976        };
1977
1978        // Inward-direction boundary check: every method the host
1979        // declares on `OdeSystem`'s M surface (and `is_autonomous`)
1980        // must round-trip through the `ParametricOdeSystem` impl.
1981        assert!(
1982            wrapper.is_autonomous(),
1983            "is_autonomous not forwarded by IcAsParametric"
1984        );
1985        assert!(
1986            wrapper.has_mass_matrix(),
1987            "has_mass_matrix not forwarded by IcAsParametric"
1988        );
1989        assert!(wrapper.is_singular_mass(), "is_singular_mass not forwarded");
1990        assert_eq!(
1991            wrapper.algebraic_indices(),
1992            vec![2],
1993            "algebraic_indices not forwarded by IcAsParametric",
1994        );
1995
1996        let mut m = vec![0.0_f64; n * n];
1997        wrapper.mass_matrix(&mut m);
1998        assert_eq!(m[0 * 3 + 0], 2.0, "mass_matrix[0,0] not forwarded");
1999        assert_eq!(m[1 * 3 + 1], 3.0, "mass_matrix[1,1] not forwarded");
2000        assert_eq!(
2001            m[2 * 3 + 2],
2002            0.0,
2003            "mass_matrix[2,2] (algebraic) not forwarded"
2004        );
2005        // Off-diagonal entries must remain zero (the host M is diagonal).
2006        assert_eq!(m[0 * 3 + 1], 0.0);
2007        assert_eq!(m[1 * 3 + 2], 0.0);
2008    }
2009}