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::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/// Wraps a [`ParametricOdeSystem`] as an [`OdeSystem`] over the augmented
284/// state `z = [y; vec(S)]` (column-major sensitivity flattening).
285///
286/// The augmented dimension is `N · (1 + N_s)`. The augmented Jacobian is
287/// `block_diag(J_y, J_y, …, J_y)` (CVODES simultaneous-corrector form).
288pub struct AugmentedSystem<S: Scalar, Sys: ParametricOdeSystem<S>> {
289 /// The underlying parameterised system.
290 pub system: Sys,
291 jy_scratch: std::cell::RefCell<Vec<S>>,
292 jp_scratch: std::cell::RefCell<Vec<S>>,
293 // Reusable FD scratch lifted out of the trait's per-call allocations.
294 // Touched only on the FD path (when has_analytical_jacobian_{y,p}
295 // returns false). Each RefCell is borrowed in isolation; no overlap.
296 fd_f0: std::cell::RefCell<Vec<S>>,
297 fd_f1: std::cell::RefCell<Vec<S>>,
298 fd_y_pert: std::cell::RefCell<Vec<S>>,
299 fd_p_pert: std::cell::RefCell<Vec<S>>,
300 // Debug-only: tracks whether the first-call analytical-vs-FD
301 // consistency check has run. See `check_jacobian_flags`.
302 #[cfg(debug_assertions)]
303 flag_check_done: std::cell::Cell<bool>,
304}
305
306impl<S: Scalar, Sys: ParametricOdeSystem<S>> AugmentedSystem<S, Sys> {
307 /// Wrap a [`ParametricOdeSystem`] into an `OdeSystem`-shaped augmented
308 /// system suitable for any [`crate::Solver`].
309 pub fn new(system: Sys) -> Self {
310 let n = system.n_states();
311 let np = system.n_params();
312 Self {
313 system,
314 jy_scratch: std::cell::RefCell::new(vec![S::ZERO; n * n]),
315 jp_scratch: std::cell::RefCell::new(vec![S::ZERO; n * np]),
316 fd_f0: std::cell::RefCell::new(vec![S::ZERO; n]),
317 fd_f1: std::cell::RefCell::new(vec![S::ZERO; n]),
318 fd_y_pert: std::cell::RefCell::new(vec![S::ZERO; n]),
319 fd_p_pert: std::cell::RefCell::new(vec![S::ZERO; np]),
320 #[cfg(debug_assertions)]
321 flag_check_done: std::cell::Cell::new(false),
322 }
323 }
324
325 /// Debug-build consistency check for the analytical-Jacobian flags.
326 ///
327 /// Catches the silent-misconfiguration class of bug where a user
328 /// implements `jacobian_y` (or `jacobian_p`) analytically but forgets
329 /// to override `has_analytical_jacobian_y` (or `_p`) — in which case
330 /// the augmented hot path bypasses the analytical override entirely
331 /// and runs FD instead. The user sees "Numra is slow on my problem"
332 /// with no diagnostic path.
333 ///
334 /// On the first `rhs` call (and only then), if a flag returns `false`,
335 /// we evaluate `system.jacobian_*` (which dispatches to the user's
336 /// override or the trait's FD default) and compare against the inline
337 /// FD result. If the user did not override the method, both calls
338 /// compute FD and agree to floating-point noise. If the user did
339 /// override with analytical, the two will typically disagree well
340 /// beyond the relative threshold — and we panic with a clear message.
341 ///
342 /// The threshold is `1e-3` relative (Frobenius). FD is only accurate
343 /// to `~sqrt(eps_mach) ≈ 1.5e-8`, so 1e-3 is loose enough to absorb
344 /// honest FD-vs-analytical disagreement at the initial state without
345 /// false positives, while still catching the "forgot the flag"
346 /// scenario where the analytical and FD results differ by O(1).
347 ///
348 /// Wrapped in `#[cfg(debug_assertions)]` — release builds pay zero.
349 #[cfg(debug_assertions)]
350 fn check_jacobian_flags(&self, t: S, y: &[S]) {
351 let n = self.system.n_states();
352 let np = self.system.n_params();
353 let threshold = S::from_f64(1e-3);
354
355 if !self.system.has_analytical_jacobian_y() {
356 let mut user_jy = vec![S::ZERO; n * n];
357 let mut fd_jy = vec![S::ZERO; n * n];
358 self.system.jacobian_y(t, y, &mut user_jy);
359 self.fd_jacobian_y_inline(t, y, &mut fd_jy);
360 if jacobians_differ_significantly(&user_jy, &fd_jy, threshold) {
361 panic!(
362 "ParametricOdeSystem implements `jacobian_y` analytically \
363 but `has_analytical_jacobian_y()` returns `false`; the \
364 analytical implementation is being ignored and `AugmentedSystem` \
365 is running finite differences instead. Override \
366 `has_analytical_jacobian_y` to return `true`. \
367 (This check is debug-build only; release builds will \
368 silently use FD.)"
369 );
370 }
371 }
372
373 if !self.system.has_analytical_jacobian_p() && np > 0 {
374 let mut user_jp = vec![S::ZERO; n * np];
375 let mut fd_jp = vec![S::ZERO; n * np];
376 self.system.jacobian_p(t, y, &mut user_jp);
377 self.fd_jacobian_p_inline(t, y, &mut fd_jp);
378 if jacobians_differ_significantly(&user_jp, &fd_jp, threshold) {
379 panic!(
380 "ParametricOdeSystem implements `jacobian_p` analytically \
381 but `has_analytical_jacobian_p()` returns `false`; the \
382 analytical implementation is being ignored and `AugmentedSystem` \
383 is running finite differences instead. Override \
384 `has_analytical_jacobian_p` to return `true`. \
385 (This check is debug-build only; release builds will \
386 silently use FD.)"
387 );
388 }
389 }
390 }
391
392 /// Forward-FD `J_y` using AugmentedSystem-local scratch buffers. Only
393 /// invoked when `system.has_analytical_jacobian_y()` is `false`.
394 fn fd_jacobian_y_inline(&self, t: S, y: &[S], jy: &mut [S]) {
395 let n = self.system.n_states();
396 let p = self.system.params();
397 let h_factor = S::EPSILON.sqrt();
398 let mut f0 = self.fd_f0.borrow_mut();
399 let mut f1 = self.fd_f1.borrow_mut();
400 let mut y_pert = self.fd_y_pert.borrow_mut();
401
402 self.system.rhs_with_params(t, y, p, &mut f0);
403 y_pert.copy_from_slice(y);
404 for j in 0..n {
405 let yj = y_pert[j];
406 let h = h_factor * (S::ONE + yj.abs());
407 y_pert[j] = yj + h;
408 self.system.rhs_with_params(t, &y_pert, p, &mut f1);
409 y_pert[j] = yj;
410 for i in 0..n {
411 jy[i * n + j] = (f1[i] - f0[i]) / h;
412 }
413 }
414 }
415
416 /// Forward-FD `J_p` using AugmentedSystem-local scratch buffers. Only
417 /// invoked when `system.has_analytical_jacobian_p()` is `false`.
418 fn fd_jacobian_p_inline(&self, t: S, y: &[S], jp: &mut [S]) {
419 let n = self.system.n_states();
420 let np = self.system.n_params();
421 let p_nominal = self.system.params();
422 let h_factor = S::EPSILON.sqrt();
423 let mut f0 = self.fd_f0.borrow_mut();
424 let mut f1 = self.fd_f1.borrow_mut();
425 let mut p_pert = self.fd_p_pert.borrow_mut();
426
427 self.system.rhs_with_params(t, y, p_nominal, &mut f0);
428 p_pert.copy_from_slice(p_nominal);
429 for k in 0..np {
430 let pk = p_pert[k];
431 let h = h_factor * (S::ONE + pk.abs());
432 p_pert[k] = pk + h;
433 self.system.rhs_with_params(t, y, &p_pert, &mut f1);
434 p_pert[k] = pk;
435 for i in 0..n {
436 jp[k * n + i] = (f1[i] - f0[i]) / h;
437 }
438 }
439 }
440
441 /// Total augmented dimension `N · (1 + N_s)`.
442 pub fn augmented_dim(&self) -> usize {
443 self.system.n_states() * (1 + self.system.n_params())
444 }
445
446 /// Build the initial augmented state `[y_0; vec(∂y_0/∂p)]`.
447 pub fn initial_augmented(&self, y0: &[S]) -> Vec<S> {
448 let n = self.system.n_states();
449 let np = self.system.n_params();
450
451 let mut z0 = Vec::with_capacity(n * (1 + np));
452 z0.extend_from_slice(y0);
453
454 let mut s0 = vec![S::ZERO; n * np];
455 self.system.initial_sensitivity(y0, &mut s0);
456 z0.extend_from_slice(&s0);
457 z0
458 }
459}
460
461/// Frobenius-norm relative comparison: returns `true` if
462/// `||a - b||_F / max(||b||_F, 1) > threshold`. Used by the debug-build
463/// analytical-vs-FD Jacobian flag check.
464#[cfg(debug_assertions)]
465fn jacobians_differ_significantly<S: Scalar>(a: &[S], b: &[S], threshold: S) -> bool {
466 debug_assert_eq!(a.len(), b.len());
467 let mut diff_sq = S::ZERO;
468 let mut b_sq = S::ZERO;
469 for i in 0..a.len() {
470 let d = a[i] - b[i];
471 diff_sq = diff_sq + d * d;
472 b_sq = b_sq + b[i] * b[i];
473 }
474 let denom = b_sq.sqrt().max(S::ONE);
475 let rel = diff_sq.sqrt() / denom;
476 rel > threshold
477}
478
479impl<S: Scalar, Sys: ParametricOdeSystem<S>> OdeSystem<S> for AugmentedSystem<S, Sys> {
480 fn dim(&self) -> usize {
481 self.augmented_dim()
482 }
483
484 /// Augmented RHS: `dy/dt = f(t,y,p)`, `dS_{:,k}/dt = J_y · S_{:,k} + J_p_{:,k}`.
485 fn rhs(&self, t: S, z: &[S], dz: &mut [S]) {
486 let n = self.system.n_states();
487 let np = self.system.n_params();
488 let y = &z[..n];
489
490 // Debug-only one-time consistency check: analytical-Jacobian flags
491 // vs the user's actual `jacobian_*` override. See
492 // `check_jacobian_flags` for rationale.
493 #[cfg(debug_assertions)]
494 if !self.flag_check_done.get() {
495 self.flag_check_done.set(true);
496 self.check_jacobian_flags(t, y);
497 }
498
499 // (a) original dynamics.
500 self.system.rhs(t, y, &mut dz[..n]);
501
502 // (b) Jacobians at the current state. Analytical when the system
503 // flags it; otherwise inline FD using local scratch buffers (lifted
504 // out of the trait defaults so the integration's hot path doesn't
505 // pay per-call allocation).
506 let mut jy = self.jy_scratch.borrow_mut();
507 let mut jp = self.jp_scratch.borrow_mut();
508 if self.system.has_analytical_jacobian_y() {
509 self.system.jacobian_y(t, y, &mut jy);
510 } else {
511 drop(jy);
512 self.fd_jacobian_y_inline(t, y, &mut self.jy_scratch.borrow_mut());
513 jy = self.jy_scratch.borrow_mut();
514 }
515 if self.system.has_analytical_jacobian_p() {
516 self.system.jacobian_p(t, y, &mut jp);
517 } else {
518 drop(jp);
519 self.fd_jacobian_p_inline(t, y, &mut self.jp_scratch.borrow_mut());
520 jp = self.jp_scratch.borrow_mut();
521 }
522
523 // (c) per-parameter sensitivity column: dS_{:,k}/dt = J_y · S_{:,k} + J_p_{:,k}.
524 // Sensitivity column-major: S_{i,k} lives at z[n + k*n + i].
525 for k in 0..np {
526 for i in 0..n {
527 let mut acc = S::ZERO;
528 // J_y row-major dot S_{:,k} (contiguous slice of z).
529 for j in 0..n {
530 acc = acc + jy[i * n + j] * z[n + k * n + j];
531 }
532 acc = acc + jp[k * n + i];
533 dz[n + k * n + i] = acc;
534 }
535 }
536 }
537
538 /// Augmented Jacobian, row-major over the full `N(1+N_s) × N(1+N_s)` block.
539 ///
540 /// Fills the `N_s + 1` diagonal blocks with `J_y`. Off-diagonal blocks
541 /// (which would carry `∂(J_y · S_{:,k})/∂y` second-derivative coupling)
542 /// are zero — the standard CVODES *simultaneous-corrector* approximation.
543 /// Implicit solvers in v1 perform dense LU on the full block; the
544 /// block-diagonal-aware fast path is a tracked follow-up.
545 fn jacobian(&self, t: S, z: &[S], jac: &mut [S]) {
546 let n = self.system.n_states();
547 let np = self.system.n_params();
548 let dim = n * (1 + np);
549 let y = &z[..n];
550
551 // Zero the entire augmented matrix first.
552 for slot in jac.iter_mut() {
553 *slot = S::ZERO;
554 }
555
556 let mut jy = self.jy_scratch.borrow_mut();
557 self.system.jacobian_y(t, y, &mut jy);
558
559 // Place J_y on each of the n_s + 1 diagonal blocks.
560 for block in 0..(1 + np) {
561 let row0 = block * n;
562 let col0 = block * n;
563 for i in 0..n {
564 for j in 0..n {
565 jac[(row0 + i) * dim + (col0 + j)] = jy[i * n + j];
566 }
567 }
568 }
569 }
570}
571
572/// Forwarding impl: a reference to a parametric system is itself a parametric
573/// system. Lets [`solve_forward_sensitivity`] accept `&Sys` without taking
574/// ownership of the user's data.
575impl<S: Scalar, T: ParametricOdeSystem<S>> ParametricOdeSystem<S> for &T {
576 fn n_states(&self) -> usize {
577 (*self).n_states()
578 }
579 fn n_params(&self) -> usize {
580 (*self).n_params()
581 }
582 fn params(&self) -> &[S] {
583 (*self).params()
584 }
585 fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]) {
586 (*self).rhs_with_params(t, y, p, dydt)
587 }
588 fn rhs(&self, t: S, y: &[S], dydt: &mut [S]) {
589 (*self).rhs(t, y, dydt)
590 }
591 fn jacobian_y(&self, t: S, y: &[S], jac: &mut [S]) {
592 (*self).jacobian_y(t, y, jac)
593 }
594 fn jacobian_p(&self, t: S, y: &[S], jp: &mut [S]) {
595 (*self).jacobian_p(t, y, jp)
596 }
597 fn initial_sensitivity(&self, y0: &[S], s0: &mut [S]) {
598 (*self).initial_sensitivity(y0, s0)
599 }
600 fn has_analytical_jacobian_y(&self) -> bool {
601 (*self).has_analytical_jacobian_y()
602 }
603 fn has_analytical_jacobian_p(&self) -> bool {
604 (*self).has_analytical_jacobian_p()
605 }
606}
607
608/// Result of an ODE forward-sensitivity solve.
609///
610/// Layout matches the conventions in the [module-level docs](self):
611/// the state `y` is row-major over time; the sensitivity is row-major over
612/// time, column-major within each per-time block.
613#[derive(Clone, Debug)]
614pub struct SensitivityResult<S: Scalar> {
615 /// Output time points (length `n_times`).
616 pub t: Vec<S>,
617 /// State trajectory: `y[i*N + j] = y_j(t_i)`, length `n_times * N`.
618 pub y: Vec<S>,
619 /// Sensitivities: per-time block of length `N · N_s`, column-major
620 /// within the block, so
621 /// `sensitivity[i * (N*N_s) + k*N + j] = ∂y_j(t_i)/∂p_k`.
622 pub sensitivity: Vec<S>,
623 /// Number of state variables `N`.
624 pub n_states: usize,
625 /// Number of parameters `N_s`.
626 pub n_params: usize,
627 /// Solver statistics from the underlying integration.
628 pub stats: SolverStats,
629 /// Was the integration successful?
630 pub success: bool,
631 /// Error description if the integration failed.
632 pub message: String,
633}
634
635impl<S: Scalar> SensitivityResult<S> {
636 /// Number of output time points.
637 pub fn len(&self) -> usize {
638 self.t.len()
639 }
640
641 /// `true` if there are no output time points.
642 pub fn is_empty(&self) -> bool {
643 self.t.is_empty()
644 }
645
646 /// State vector at output index `i` (length `N`).
647 pub fn y_at(&self, i: usize) -> &[S] {
648 let n = self.n_states;
649 let off = i * n;
650 &self.y[off..off + n]
651 }
652
653 /// Full sensitivity block at output index `i`, column-major
654 /// (length `N · N_s`). `sensitivity_at(i)[k*N + j] = ∂y_j(t_i)/∂p_k`.
655 pub fn sensitivity_at(&self, i: usize) -> &[S] {
656 let block = self.n_states * self.n_params;
657 let off = i * block;
658 &self.sensitivity[off..off + block]
659 }
660
661 /// Contiguous slice `S_{:,k}(t_i) = ∂y/∂p_k` at output index `i`,
662 /// length `N`. Free slice thanks to the column-major flattening.
663 pub fn sensitivity_for_param(&self, i: usize, k: usize) -> &[S] {
664 let n = self.n_states;
665 let block = n * self.n_params;
666 let off = i * block + k * n;
667 &self.sensitivity[off..off + n]
668 }
669
670 /// Single component `∂y_state(t_i)/∂p_param`.
671 pub fn dyi_dpj(&self, i: usize, state: usize, param: usize) -> S {
672 let n = self.n_states;
673 let block = n * self.n_params;
674 self.sensitivity[i * block + param * n + state]
675 }
676
677 /// State at the final output time. Returns an empty slice if the result
678 /// has no time points.
679 pub fn final_state(&self) -> &[S] {
680 if self.is_empty() {
681 &[]
682 } else {
683 self.y_at(self.len() - 1)
684 }
685 }
686
687 /// Sensitivity block at the final output time, column-major (length
688 /// `N · N_s`). Empty slice when the result has no time points.
689 pub fn final_sensitivity(&self) -> &[S] {
690 if self.is_empty() {
691 &[]
692 } else {
693 self.sensitivity_at(self.len() - 1)
694 }
695 }
696
697 /// Normalised (logarithmic) sensitivity at output index `i`:
698 /// `(p_k / y_j) · ∂y_j/∂p_k`, column-major in the same shape as
699 /// [`Self::sensitivity_at`]. Components where `|y_j|` falls below
700 /// `1e-15 · max(1, |y|_∞)` are reported as zero to avoid blow-up.
701 pub fn normalized_sensitivity_at(&self, i: usize, p_nominal: &[S]) -> Vec<S> {
702 assert_eq!(
703 p_nominal.len(),
704 self.n_params,
705 "p_nominal length {} does not match n_params {}",
706 p_nominal.len(),
707 self.n_params,
708 );
709 let n = self.n_states;
710 let np = self.n_params;
711 let y = self.y_at(i);
712 let s = self.sensitivity_at(i);
713
714 let mut y_max = S::ZERO;
715 for &v in y {
716 let av = v.abs();
717 if av > y_max {
718 y_max = av;
719 }
720 }
721 let zero_threshold = S::from_f64(1e-15) * (S::ONE + y_max);
722
723 let mut out = vec![S::ZERO; n * np];
724 for k in 0..np {
725 for j in 0..n {
726 let yj = y[j];
727 if yj.abs() <= zero_threshold {
728 out[k * n + j] = S::ZERO;
729 } else {
730 out[k * n + j] = (p_nominal[k] / yj) * s[k * n + j];
731 }
732 }
733 }
734 out
735 }
736}
737
738/// Closure-shaped adapter for [`solve_forward_sensitivity_with`].
739///
740/// Wraps a closure of shape `Fn(t, y, p, dydt)` and a fixed parameter vector
741/// into a [`ParametricOdeSystem`] using FD-default Jacobians. Public so that
742/// callers writing a one-shot harness can construct it explicitly when they
743/// need to share it across multiple solver invocations.
744pub struct ClosureSystem<S: Scalar, F> {
745 rhs: F,
746 params: Vec<S>,
747 n_states: usize,
748}
749
750impl<S: Scalar, F> ClosureSystem<S, F>
751where
752 F: Fn(S, &[S], &[S], &mut [S]),
753{
754 /// Build a closure-backed [`ParametricOdeSystem`].
755 pub fn new(rhs: F, params: Vec<S>, n_states: usize) -> Self {
756 Self {
757 rhs,
758 params,
759 n_states,
760 }
761 }
762}
763
764impl<S: Scalar, F> ParametricOdeSystem<S> for ClosureSystem<S, F>
765where
766 F: Fn(S, &[S], &[S], &mut [S]),
767{
768 fn n_states(&self) -> usize {
769 self.n_states
770 }
771 fn n_params(&self) -> usize {
772 self.params.len()
773 }
774 fn params(&self) -> &[S] {
775 &self.params
776 }
777 fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]) {
778 (self.rhs)(t, y, p, dydt)
779 }
780}
781
782/// Compute the forward sensitivity matrix `S(t) = ∂y(t) / ∂p` of an ODE
783/// solution using any [`Solver`].
784///
785/// Builds the augmented state `z = [y; vec(S)]` (column-major sensitivity
786/// flattening; see [module-level docs](self)), drives `Sol` over the
787/// augmented system, and splits the trajectory back into `y` and
788/// sensitivity blocks at every output time the solver produced.
789///
790/// # Errors
791///
792/// Returns whatever error `Sol::solve` raises. If the underlying integration
793/// returns `success = false` (e.g. step-size collapse or max-steps reached),
794/// the returned [`SensitivityResult`] propagates that with empty trajectory
795/// vectors and the solver's diagnostic message.
796///
797/// # Example
798///
799/// ```
800/// use numra_ode::sensitivity::{solve_forward_sensitivity, ParametricOdeSystem};
801/// use numra_ode::{DoPri5, SolverOptions};
802///
803/// // dy/dt = -k * y, y(0) = 1, k = 0.5.
804/// // Analytical: y(t) = exp(-k t), ∂y/∂k = -t exp(-k t).
805/// struct Decay { k: f64 }
806/// impl ParametricOdeSystem<f64> for Decay {
807/// fn n_states(&self) -> usize { 1 }
808/// fn n_params(&self) -> usize { 1 }
809/// fn params(&self) -> &[f64] { std::slice::from_ref(&self.k) }
810/// fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
811/// dy[0] = -p[0] * y[0];
812/// }
813/// fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) { jy[0] = -self.k; }
814/// fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) { jp[0] = -y[0]; }
815/// fn has_analytical_jacobian_y(&self) -> bool { true }
816/// fn has_analytical_jacobian_p(&self) -> bool { true }
817/// }
818///
819/// let r = solve_forward_sensitivity::<DoPri5, _, _>(
820/// &Decay { k: 0.5 },
821/// 0.0, 2.0,
822/// &[1.0],
823/// &SolverOptions::default().rtol(1e-8).atol(1e-10),
824/// ).unwrap();
825///
826/// let last = r.len() - 1;
827/// let dy_dk = r.dyi_dpj(last, 0, 0);
828/// let analytical = -r.t[last] * (-0.5 * r.t[last]).exp();
829/// assert!((dy_dk - analytical).abs() < 1e-5);
830/// ```
831pub fn solve_forward_sensitivity<Sol, S, Sys>(
832 system: &Sys,
833 t0: S,
834 tf: S,
835 y0: &[S],
836 options: &SolverOptions<S>,
837) -> Result<SensitivityResult<S>, SolverError>
838where
839 Sol: Solver<S>,
840 S: Scalar,
841 Sys: ParametricOdeSystem<S>,
842{
843 let n = system.n_states();
844 let np = system.n_params();
845 assert_eq!(
846 y0.len(),
847 n,
848 "solve_forward_sensitivity: y0.len() = {} but n_states = {}",
849 y0.len(),
850 n,
851 );
852 assert_eq!(
853 system.params().len(),
854 np,
855 "solve_forward_sensitivity: params().len() = {} but n_params = {}",
856 system.params().len(),
857 np,
858 );
859
860 let aug = AugmentedSystem::new(system);
861 let z0 = aug.initial_augmented(y0);
862 let aug_dim = aug.augmented_dim();
863
864 let aug_result = Sol::solve(&aug, t0, tf, &z0, options)?;
865
866 if !aug_result.success {
867 return Ok(SensitivityResult {
868 t: Vec::new(),
869 y: Vec::new(),
870 sensitivity: Vec::new(),
871 n_states: n,
872 n_params: np,
873 stats: aug_result.stats,
874 success: false,
875 message: aug_result.message,
876 });
877 }
878
879 let n_times = aug_result.len();
880 let mut t = Vec::with_capacity(n_times);
881 let mut y = Vec::with_capacity(n_times * n);
882 let mut sensitivity = Vec::with_capacity(n_times * n * np);
883
884 for i in 0..n_times {
885 t.push(aug_result.t[i]);
886 let block_start = i * aug_dim;
887 let state_block = &aug_result.y[block_start..block_start + n];
888 let sens_block = &aug_result.y[block_start + n..block_start + aug_dim];
889 y.extend_from_slice(state_block);
890 sensitivity.extend_from_slice(sens_block);
891 }
892
893 Ok(SensitivityResult {
894 t,
895 y,
896 sensitivity,
897 n_states: n,
898 n_params: np,
899 stats: aug_result.stats,
900 success: true,
901 message: String::new(),
902 })
903}
904
905/// Closure-shaped convenience wrapper around [`solve_forward_sensitivity`].
906///
907/// Wraps `rhs(t, y, p, dydt)` plus a parameter vector into an internal
908/// [`ClosureSystem`] (FD-default Jacobians) and forwards. Suitable for
909/// one-shot analyses, tests, and REPL-style scripts. For production usage
910/// — particularly stiff problems where analytical Jacobians materially
911/// change runtime — implement [`ParametricOdeSystem`] directly and call
912/// [`solve_forward_sensitivity`].
913///
914/// # Example
915///
916/// ```
917/// use numra_ode::{DoPri5, SolverOptions};
918/// use numra_ode::sensitivity::solve_forward_sensitivity_with;
919///
920/// // dy/dt = -k * y, y(0) = 1, k = 0.5.
921/// let r = solve_forward_sensitivity_with::<DoPri5, f64, _>(
922/// |_t, y, p, dy| { dy[0] = -p[0] * y[0]; },
923/// &[1.0],
924/// &[0.5],
925/// 0.0, 2.0,
926/// &SolverOptions::default().rtol(1e-8).atol(1e-10),
927/// ).unwrap();
928///
929/// let last = r.len() - 1;
930/// let dy_dk = r.dyi_dpj(last, 0, 0);
931/// let analytical: f64 = -r.t[last] * (-0.5 * r.t[last]).exp();
932/// assert!((dy_dk - analytical).abs() < 1e-4);
933/// ```
934pub fn solve_forward_sensitivity_with<Sol, S, F>(
935 rhs: F,
936 y0: &[S],
937 params: &[S],
938 t0: S,
939 tf: S,
940 options: &SolverOptions<S>,
941) -> Result<SensitivityResult<S>, SolverError>
942where
943 Sol: Solver<S>,
944 S: Scalar,
945 F: Fn(S, &[S], &[S], &mut [S]),
946{
947 let system = ClosureSystem {
948 rhs,
949 params: params.to_vec(),
950 n_states: y0.len(),
951 };
952 solve_forward_sensitivity::<Sol, S, _>(&system, t0, tf, y0, options)
953}
954
955#[cfg(test)]
956mod tests {
957 use super::*;
958
959 /// Minimal scalar decay system used to drive the trait through
960 /// `OdeSystem::rhs` and `OdeSystem::jacobian` without involving a solver.
961 struct ExpDecay {
962 k: f64,
963 }
964
965 impl ParametricOdeSystem<f64> for ExpDecay {
966 fn n_states(&self) -> usize {
967 1
968 }
969 fn n_params(&self) -> usize {
970 1
971 }
972 fn params(&self) -> &[f64] {
973 std::slice::from_ref(&self.k)
974 }
975 fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
976 dy[0] = -p[0] * y[0];
977 }
978 fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) {
979 jy[0] = -self.k;
980 }
981 fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) {
982 jp[0] = -y[0];
983 }
984 fn has_analytical_jacobian_y(&self) -> bool {
985 true
986 }
987 fn has_analytical_jacobian_p(&self) -> bool {
988 true
989 }
990 }
991
992 /// 2-state, 2-parameter linear system (mirrors the module-doc example).
993 /// Lets us exercise the asymmetric J_y / J_p layouts in code.
994 struct Lin2 {
995 p: [f64; 2],
996 }
997
998 impl ParametricOdeSystem<f64> for Lin2 {
999 fn n_states(&self) -> usize {
1000 2
1001 }
1002 fn n_params(&self) -> usize {
1003 2
1004 }
1005 fn params(&self) -> &[f64] {
1006 &self.p
1007 }
1008 fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
1009 dy[0] = -p[0] * y[0] + p[1] * y[1];
1010 dy[1] = -p[1] * y[1];
1011 }
1012 fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) {
1013 jy[0 * 2 + 0] = -self.p[0];
1014 jy[0 * 2 + 1] = self.p[1];
1015 jy[1 * 2 + 0] = 0.0;
1016 jy[1 * 2 + 1] = -self.p[1];
1017 }
1018 fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) {
1019 jp[0 * 2 + 0] = -y[0];
1020 jp[0 * 2 + 1] = 0.0;
1021 jp[1 * 2 + 0] = y[1];
1022 jp[1 * 2 + 1] = -y[1];
1023 }
1024 fn has_analytical_jacobian_y(&self) -> bool {
1025 true
1026 }
1027 fn has_analytical_jacobian_p(&self) -> bool {
1028 true
1029 }
1030 }
1031
1032 #[test]
1033 fn augmented_dim_and_initial_state() {
1034 let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
1035 assert_eq!(aug.augmented_dim(), 2 * (1 + 2));
1036 let z0 = aug.initial_augmented(&[1.0, 2.0]);
1037 assert_eq!(z0.len(), 6);
1038 assert_eq!(&z0[..2], &[1.0, 2.0]);
1039 for s in &z0[2..] {
1040 assert_eq!(*s, 0.0);
1041 }
1042 }
1043
1044 #[test]
1045 fn augmented_rhs_matches_hand_derivation() {
1046 // System: dy/dt = -k y, S_{0,0} = ∂y/∂k.
1047 // dS/dt = J_y * S + J_p = -k * S - y.
1048 let aug = AugmentedSystem::new(ExpDecay { k: 0.5 });
1049 let z = vec![2.0, 1.0]; // y = 2.0, S = 1.0
1050 let mut dz = vec![0.0; 2];
1051 aug.rhs(0.0, &z, &mut dz);
1052
1053 // dy/dt = -0.5 * 2.0 = -1.0
1054 assert!((dz[0] + 1.0).abs() < 1e-12);
1055 // dS/dt = -0.5 * 1.0 + (-2.0) = -2.5
1056 assert!((dz[1] + 2.5).abs() < 1e-12);
1057 }
1058
1059 #[test]
1060 fn augmented_rhs_matches_hand_derivation_two_param() {
1061 // Lin2 system at y = (1, 2), p = (1, 0.5), S = identity (initial sensitivity ID).
1062 let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
1063
1064 // z = [y_0, y_1, S_{0,0}, S_{1,0}, S_{0,1}, S_{1,1}]
1065 // [ 1, 2, 1, 0, 0, 1 ]
1066 let z = vec![1.0, 2.0, 1.0, 0.0, 0.0, 1.0];
1067 let mut dz = vec![0.0; 6];
1068 aug.rhs(0.0, &z, &mut dz);
1069
1070 // dy/dt
1071 // dy_0/dt = -1*1 + 0.5*2 = 0.0
1072 // dy_1/dt = -0.5*2 = -1.0
1073 assert!((dz[0] - 0.0).abs() < 1e-12);
1074 assert!((dz[1] + 1.0).abs() < 1e-12);
1075
1076 // dS_{:,0}/dt = J_y · [1; 0] + J_p_{:,0}
1077 // J_y · [1; 0] = [-1; 0]; J_p_{:,0} = [-1; 0]; total = [-2; 0]
1078 assert!((dz[2] + 2.0).abs() < 1e-12);
1079 assert!((dz[3] - 0.0).abs() < 1e-12);
1080
1081 // dS_{:,1}/dt = J_y · [0; 1] + J_p_{:,1}
1082 // J_y · [0; 1] = [0.5; -0.5]; J_p_{:,1} = [2; -2]; total = [2.5; -2.5]
1083 assert!((dz[4] - 2.5).abs() < 1e-12);
1084 assert!((dz[5] + 2.5).abs() < 1e-12);
1085 }
1086
1087 #[test]
1088 fn augmented_jacobian_is_block_diagonal_with_jy() {
1089 let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
1090 let z = vec![1.0, 2.0, 1.0, 0.0, 0.0, 1.0];
1091 let dim = aug.augmented_dim();
1092 let mut jac = vec![0.0; dim * dim];
1093 aug.jacobian(0.0, &z, &mut jac);
1094
1095 // Expected J_y at p=(1, 0.5):
1096 // [-1.0 0.5]
1097 // [ 0.0 -0.5]
1098 let jy = [-1.0_f64, 0.5, 0.0, -0.5];
1099
1100 for block in 0..3 {
1101 let row0 = block * 2;
1102 let col0 = block * 2;
1103 for i in 0..2 {
1104 for j in 0..2 {
1105 let got = jac[(row0 + i) * dim + (col0 + j)];
1106 let want = jy[i * 2 + j];
1107 assert!(
1108 (got - want).abs() < 1e-12,
1109 "diag block {block}, ({i},{j}): got {got}, want {want}"
1110 );
1111 }
1112 }
1113 }
1114
1115 // Off-diagonal blocks must be zero.
1116 for br in 0..3 {
1117 for bc in 0..3 {
1118 if br == bc {
1119 continue;
1120 }
1121 for i in 0..2 {
1122 for j in 0..2 {
1123 let r = br * 2 + i;
1124 let c = bc * 2 + j;
1125 assert_eq!(jac[r * dim + c], 0.0, "off-diag ({br},{bc}) ({i},{j})");
1126 }
1127 }
1128 }
1129 }
1130 }
1131
1132 #[test]
1133 fn fd_default_jacobian_y_matches_analytical() {
1134 // Re-implement Lin2 without analytical J_y so the FD default is exercised.
1135 struct Lin2Fd {
1136 p: [f64; 2],
1137 }
1138 impl ParametricOdeSystem<f64> for Lin2Fd {
1139 fn n_states(&self) -> usize {
1140 2
1141 }
1142 fn n_params(&self) -> usize {
1143 2
1144 }
1145 fn params(&self) -> &[f64] {
1146 &self.p
1147 }
1148 fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
1149 dy[0] = -p[0] * y[0] + p[1] * y[1];
1150 dy[1] = -p[1] * y[1];
1151 }
1152 // No jacobian_y / jacobian_p override — FD defaults are used.
1153 }
1154
1155 let analytical = Lin2 { p: [1.0, 0.5] };
1156 let fd = Lin2Fd { p: [1.0, 0.5] };
1157
1158 let y = [1.0, 2.0];
1159 let mut a_jy = [0.0; 4];
1160 let mut f_jy = [0.0; 4];
1161 analytical.jacobian_y(0.0, &y, &mut a_jy);
1162 fd.jacobian_y(0.0, &y, &mut f_jy);
1163 for k in 0..4 {
1164 assert!(
1165 (a_jy[k] - f_jy[k]).abs() < 1e-7,
1166 "j_y[{k}] analytical={} fd={}",
1167 a_jy[k],
1168 f_jy[k]
1169 );
1170 }
1171
1172 let mut a_jp = [0.0; 4];
1173 let mut f_jp = [0.0; 4];
1174 analytical.jacobian_p(0.0, &y, &mut a_jp);
1175 fd.jacobian_p(0.0, &y, &mut f_jp);
1176 for k in 0..4 {
1177 assert!(
1178 (a_jp[k] - f_jp[k]).abs() < 1e-7,
1179 "j_p[{k}] analytical={} fd={}",
1180 a_jp[k],
1181 f_jp[k]
1182 );
1183 }
1184 }
1185
1186 #[test]
1187 fn sensitivity_result_accessors() {
1188 // Hand-construct a result for a 2-state, 3-param trajectory at 2 time points.
1189 // sensitivity layout per time block (column-major, length N*N_s = 6):
1190 // [ S_{0,0} S_{1,0} | S_{0,1} S_{1,1} | S_{0,2} S_{1,2} ]
1191 let res = SensitivityResult {
1192 t: vec![0.0, 1.0],
1193 y: vec![1.0, 2.0, 0.5, 1.5],
1194 sensitivity: vec![
1195 // t = 0
1196 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, // t = 1
1197 1.1, 1.2, 1.3, 1.4, 1.5, 1.6,
1198 ],
1199 n_states: 2,
1200 n_params: 3,
1201 stats: SolverStats::new(),
1202 success: true,
1203 message: String::new(),
1204 };
1205
1206 assert_eq!(res.len(), 2);
1207 assert_eq!(res.y_at(1), &[0.5, 1.5]);
1208 assert_eq!(res.sensitivity_at(0).len(), 6);
1209 assert_eq!(res.sensitivity_for_param(1, 2), &[1.5, 1.6]);
1210 assert!((res.dyi_dpj(0, 1, 0) - 0.2).abs() < 1e-12);
1211 assert!((res.dyi_dpj(1, 0, 2) - 1.5).abs() < 1e-12);
1212 assert_eq!(res.final_state(), &[0.5, 1.5]);
1213 assert_eq!(res.final_sensitivity().len(), 6);
1214 }
1215
1216 #[test]
1217 fn normalized_sensitivity_handles_zero_state_safely() {
1218 let res = SensitivityResult {
1219 t: vec![0.0],
1220 y: vec![1.0, 0.0],
1221 sensitivity: vec![0.5, 0.7],
1222 n_states: 2,
1223 n_params: 1,
1224 stats: SolverStats::new(),
1225 success: true,
1226 message: String::new(),
1227 };
1228 let norm = res.normalized_sensitivity_at(0, &[2.0]);
1229 // (p / y_0) · S_{0,0} = (2 / 1) · 0.5 = 1.0
1230 assert!((norm[0] - 1.0).abs() < 1e-12);
1231 // y_1 = 0 → reported as zero rather than ∞.
1232 assert_eq!(norm[1], 0.0);
1233 }
1234
1235 // ---------------------------------------------------------------
1236 // Solve entry-point sanity tests.
1237 //
1238 // The full solver-coverage matrix lives in
1239 // numra-ode/tests/sensitivity_regression.rs (commit 3).
1240 // ---------------------------------------------------------------
1241
1242 use crate::DoPri5;
1243 use crate::SolverOptions;
1244
1245 #[test]
1246 fn solve_forward_sensitivity_with_dopri5_matches_analytical_decay() {
1247 // dy/dt = -k y, y(0) = 1, k = 0.5.
1248 // Analytical: ∂y/∂k(t) = -t · exp(-k t).
1249 let result = solve_forward_sensitivity_with::<DoPri5, f64, _>(
1250 |_t, y, p, dy| {
1251 dy[0] = -p[0] * y[0];
1252 },
1253 &[1.0],
1254 &[0.5],
1255 0.0,
1256 2.0,
1257 &SolverOptions::default().rtol(1e-8).atol(1e-10),
1258 )
1259 .expect("solve failed");
1260
1261 assert!(result.success);
1262 assert_eq!(result.n_states, 1);
1263 assert_eq!(result.n_params, 1);
1264 let last = result.len() - 1;
1265 let dy_dk = result.dyi_dpj(last, 0, 0);
1266 let t_last = result.t[last];
1267 let analytical = -t_last * (-0.5 * t_last).exp();
1268 assert!(
1269 (dy_dk - analytical).abs() < 1e-5,
1270 "computed {dy_dk}, analytical {analytical}, |err| = {}",
1271 (dy_dk - analytical).abs()
1272 );
1273 }
1274
1275 #[test]
1276 fn solve_trait_form_matches_closure_form() {
1277 let trait_result = solve_forward_sensitivity::<DoPri5, _, _>(
1278 &ExpDecay { k: 0.5 },
1279 0.0,
1280 2.0,
1281 &[1.0],
1282 &SolverOptions::default().rtol(1e-9).atol(1e-12),
1283 )
1284 .expect("trait solve failed");
1285
1286 let closure_result = solve_forward_sensitivity_with::<DoPri5, f64, _>(
1287 |_t, y, p, dy| {
1288 dy[0] = -p[0] * y[0];
1289 },
1290 &[1.0],
1291 &[0.5],
1292 0.0,
1293 2.0,
1294 &SolverOptions::default().rtol(1e-9).atol(1e-12),
1295 )
1296 .expect("closure solve failed");
1297
1298 // Trait form has analytical Jacobians; closure form falls back to FD.
1299 // Both must converge to the same answer at the final time, modulo FD
1300 // noise.
1301 let n_t = trait_result.len();
1302 let n_c = closure_result.len();
1303 let trait_final = trait_result.dyi_dpj(n_t - 1, 0, 0);
1304 let closure_final = closure_result.dyi_dpj(n_c - 1, 0, 0);
1305 assert!(
1306 (trait_final - closure_final).abs() < 1e-5,
1307 "trait {trait_final} vs closure {closure_final}, |err| = {}",
1308 (trait_final - closure_final).abs()
1309 );
1310 }
1311
1312 #[test]
1313 fn solve_propagates_solver_failure_as_unsuccessful_result() {
1314 // Forcing max_steps = 1 against a multi-step problem makes the
1315 // underlying solver return success = false (or an error). We assert
1316 // that solve_forward_sensitivity surfaces this without panicking.
1317 let opts = SolverOptions::default().rtol(1e-9).atol(1e-12).max_steps(1);
1318
1319 let outcome = solve_forward_sensitivity_with::<DoPri5, f64, _>(
1320 |_t, y, p, dy| {
1321 dy[0] = -p[0] * y[0];
1322 },
1323 &[1.0],
1324 &[0.5],
1325 0.0,
1326 10.0,
1327 &opts,
1328 );
1329
1330 match outcome {
1331 Ok(r) => {
1332 assert!(!r.success, "expected unsuccessful result, got success");
1333 assert!(r.t.is_empty());
1334 assert!(r.y.is_empty());
1335 assert!(r.sensitivity.is_empty());
1336 }
1337 Err(_) => {
1338 // Some solvers return Err directly for max_steps; that's also acceptable.
1339 }
1340 }
1341 }
1342}