use numra_core::Scalar;
use crate::error::SolverError;
use crate::problem::OdeSystem;
use crate::solver::{Solver, SolverOptions, SolverStats};
pub trait ParametricOdeSystem<S: Scalar> {
fn n_states(&self) -> usize;
fn n_params(&self) -> usize;
fn params(&self) -> &[S];
fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]);
fn rhs(&self, t: S, y: &[S], dydt: &mut [S]) {
self.rhs_with_params(t, y, self.params(), dydt);
}
fn jacobian_y(&self, t: S, y: &[S], jac: &mut [S]) {
let n = self.n_states();
let p = self.params();
let h_factor = S::EPSILON.sqrt();
let mut f0 = vec![S::ZERO; n];
let mut f1 = vec![S::ZERO; n];
let mut y_pert = y.to_vec();
self.rhs_with_params(t, y, p, &mut f0);
for j in 0..n {
let yj = y_pert[j];
let h = h_factor * (S::ONE + yj.abs());
y_pert[j] = yj + h;
self.rhs_with_params(t, &y_pert, p, &mut f1);
y_pert[j] = yj;
for i in 0..n {
jac[i * n + j] = (f1[i] - f0[i]) / h;
}
}
}
fn jacobian_p(&self, t: S, y: &[S], jp: &mut [S]) {
let n = self.n_states();
let np = self.n_params();
let p_nominal = self.params();
let h_factor = S::EPSILON.sqrt();
let mut f0 = vec![S::ZERO; n];
let mut f1 = vec![S::ZERO; n];
let mut p_pert = p_nominal.to_vec();
self.rhs_with_params(t, y, p_nominal, &mut f0);
for k in 0..np {
let pk = p_pert[k];
let h = h_factor * (S::ONE + pk.abs());
p_pert[k] = pk + h;
self.rhs_with_params(t, y, &p_pert, &mut f1);
p_pert[k] = pk;
for i in 0..n {
jp[k * n + i] = (f1[i] - f0[i]) / h;
}
}
}
fn initial_sensitivity(&self, _y0: &[S], s0: &mut [S]) {
for s in s0.iter_mut() {
*s = S::ZERO;
}
}
fn has_analytical_jacobian_y(&self) -> bool {
false
}
fn has_analytical_jacobian_p(&self) -> bool {
false
}
}
pub struct AugmentedSystem<S: Scalar, Sys: ParametricOdeSystem<S>> {
pub system: Sys,
jy_scratch: std::cell::RefCell<Vec<S>>,
jp_scratch: std::cell::RefCell<Vec<S>>,
fd_f0: std::cell::RefCell<Vec<S>>,
fd_f1: std::cell::RefCell<Vec<S>>,
fd_y_pert: std::cell::RefCell<Vec<S>>,
fd_p_pert: std::cell::RefCell<Vec<S>>,
#[cfg(debug_assertions)]
flag_check_done: std::cell::Cell<bool>,
}
impl<S: Scalar, Sys: ParametricOdeSystem<S>> AugmentedSystem<S, Sys> {
pub fn new(system: Sys) -> Self {
let n = system.n_states();
let np = system.n_params();
Self {
system,
jy_scratch: std::cell::RefCell::new(vec![S::ZERO; n * n]),
jp_scratch: std::cell::RefCell::new(vec![S::ZERO; n * np]),
fd_f0: std::cell::RefCell::new(vec![S::ZERO; n]),
fd_f1: std::cell::RefCell::new(vec![S::ZERO; n]),
fd_y_pert: std::cell::RefCell::new(vec![S::ZERO; n]),
fd_p_pert: std::cell::RefCell::new(vec![S::ZERO; np]),
#[cfg(debug_assertions)]
flag_check_done: std::cell::Cell::new(false),
}
}
#[cfg(debug_assertions)]
fn check_jacobian_flags(&self, t: S, y: &[S]) {
let n = self.system.n_states();
let np = self.system.n_params();
let threshold = S::from_f64(1e-3);
if !self.system.has_analytical_jacobian_y() {
let mut user_jy = vec![S::ZERO; n * n];
let mut fd_jy = vec![S::ZERO; n * n];
self.system.jacobian_y(t, y, &mut user_jy);
self.fd_jacobian_y_inline(t, y, &mut fd_jy);
if jacobians_differ_significantly(&user_jy, &fd_jy, threshold) {
panic!(
"ParametricOdeSystem implements `jacobian_y` analytically \
but `has_analytical_jacobian_y()` returns `false`; the \
analytical implementation is being ignored and `AugmentedSystem` \
is running finite differences instead. Override \
`has_analytical_jacobian_y` to return `true`. \
(This check is debug-build only; release builds will \
silently use FD.)"
);
}
}
if !self.system.has_analytical_jacobian_p() && np > 0 {
let mut user_jp = vec![S::ZERO; n * np];
let mut fd_jp = vec![S::ZERO; n * np];
self.system.jacobian_p(t, y, &mut user_jp);
self.fd_jacobian_p_inline(t, y, &mut fd_jp);
if jacobians_differ_significantly(&user_jp, &fd_jp, threshold) {
panic!(
"ParametricOdeSystem implements `jacobian_p` analytically \
but `has_analytical_jacobian_p()` returns `false`; the \
analytical implementation is being ignored and `AugmentedSystem` \
is running finite differences instead. Override \
`has_analytical_jacobian_p` to return `true`. \
(This check is debug-build only; release builds will \
silently use FD.)"
);
}
}
}
fn fd_jacobian_y_inline(&self, t: S, y: &[S], jy: &mut [S]) {
let n = self.system.n_states();
let p = self.system.params();
let h_factor = S::EPSILON.sqrt();
let mut f0 = self.fd_f0.borrow_mut();
let mut f1 = self.fd_f1.borrow_mut();
let mut y_pert = self.fd_y_pert.borrow_mut();
self.system.rhs_with_params(t, y, p, &mut f0);
y_pert.copy_from_slice(y);
for j in 0..n {
let yj = y_pert[j];
let h = h_factor * (S::ONE + yj.abs());
y_pert[j] = yj + h;
self.system.rhs_with_params(t, &y_pert, p, &mut f1);
y_pert[j] = yj;
for i in 0..n {
jy[i * n + j] = (f1[i] - f0[i]) / h;
}
}
}
fn fd_jacobian_p_inline(&self, t: S, y: &[S], jp: &mut [S]) {
let n = self.system.n_states();
let np = self.system.n_params();
let p_nominal = self.system.params();
let h_factor = S::EPSILON.sqrt();
let mut f0 = self.fd_f0.borrow_mut();
let mut f1 = self.fd_f1.borrow_mut();
let mut p_pert = self.fd_p_pert.borrow_mut();
self.system.rhs_with_params(t, y, p_nominal, &mut f0);
p_pert.copy_from_slice(p_nominal);
for k in 0..np {
let pk = p_pert[k];
let h = h_factor * (S::ONE + pk.abs());
p_pert[k] = pk + h;
self.system.rhs_with_params(t, y, &p_pert, &mut f1);
p_pert[k] = pk;
for i in 0..n {
jp[k * n + i] = (f1[i] - f0[i]) / h;
}
}
}
pub fn augmented_dim(&self) -> usize {
self.system.n_states() * (1 + self.system.n_params())
}
pub fn initial_augmented(&self, y0: &[S]) -> Vec<S> {
let n = self.system.n_states();
let np = self.system.n_params();
let mut z0 = Vec::with_capacity(n * (1 + np));
z0.extend_from_slice(y0);
let mut s0 = vec![S::ZERO; n * np];
self.system.initial_sensitivity(y0, &mut s0);
z0.extend_from_slice(&s0);
z0
}
}
#[cfg(debug_assertions)]
fn jacobians_differ_significantly<S: Scalar>(a: &[S], b: &[S], threshold: S) -> bool {
debug_assert_eq!(a.len(), b.len());
let mut diff_sq = S::ZERO;
let mut b_sq = S::ZERO;
for i in 0..a.len() {
let d = a[i] - b[i];
diff_sq = diff_sq + d * d;
b_sq = b_sq + b[i] * b[i];
}
let denom = b_sq.sqrt().max(S::ONE);
let rel = diff_sq.sqrt() / denom;
rel > threshold
}
impl<S: Scalar, Sys: ParametricOdeSystem<S>> OdeSystem<S> for AugmentedSystem<S, Sys> {
fn dim(&self) -> usize {
self.augmented_dim()
}
fn rhs(&self, t: S, z: &[S], dz: &mut [S]) {
let n = self.system.n_states();
let np = self.system.n_params();
let y = &z[..n];
#[cfg(debug_assertions)]
if !self.flag_check_done.get() {
self.flag_check_done.set(true);
self.check_jacobian_flags(t, y);
}
self.system.rhs(t, y, &mut dz[..n]);
let mut jy = self.jy_scratch.borrow_mut();
let mut jp = self.jp_scratch.borrow_mut();
if self.system.has_analytical_jacobian_y() {
self.system.jacobian_y(t, y, &mut jy);
} else {
drop(jy);
self.fd_jacobian_y_inline(t, y, &mut self.jy_scratch.borrow_mut());
jy = self.jy_scratch.borrow_mut();
}
if self.system.has_analytical_jacobian_p() {
self.system.jacobian_p(t, y, &mut jp);
} else {
drop(jp);
self.fd_jacobian_p_inline(t, y, &mut self.jp_scratch.borrow_mut());
jp = self.jp_scratch.borrow_mut();
}
for k in 0..np {
for i in 0..n {
let mut acc = S::ZERO;
for j in 0..n {
acc = acc + jy[i * n + j] * z[n + k * n + j];
}
acc = acc + jp[k * n + i];
dz[n + k * n + i] = acc;
}
}
}
fn jacobian(&self, t: S, z: &[S], jac: &mut [S]) {
let n = self.system.n_states();
let np = self.system.n_params();
let dim = n * (1 + np);
let y = &z[..n];
for slot in jac.iter_mut() {
*slot = S::ZERO;
}
let mut jy = self.jy_scratch.borrow_mut();
self.system.jacobian_y(t, y, &mut jy);
for block in 0..(1 + np) {
let row0 = block * n;
let col0 = block * n;
for i in 0..n {
for j in 0..n {
jac[(row0 + i) * dim + (col0 + j)] = jy[i * n + j];
}
}
}
}
}
impl<S: Scalar, T: ParametricOdeSystem<S>> ParametricOdeSystem<S> for &T {
fn n_states(&self) -> usize {
(*self).n_states()
}
fn n_params(&self) -> usize {
(*self).n_params()
}
fn params(&self) -> &[S] {
(*self).params()
}
fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]) {
(*self).rhs_with_params(t, y, p, dydt)
}
fn rhs(&self, t: S, y: &[S], dydt: &mut [S]) {
(*self).rhs(t, y, dydt)
}
fn jacobian_y(&self, t: S, y: &[S], jac: &mut [S]) {
(*self).jacobian_y(t, y, jac)
}
fn jacobian_p(&self, t: S, y: &[S], jp: &mut [S]) {
(*self).jacobian_p(t, y, jp)
}
fn initial_sensitivity(&self, y0: &[S], s0: &mut [S]) {
(*self).initial_sensitivity(y0, s0)
}
fn has_analytical_jacobian_y(&self) -> bool {
(*self).has_analytical_jacobian_y()
}
fn has_analytical_jacobian_p(&self) -> bool {
(*self).has_analytical_jacobian_p()
}
}
#[derive(Clone, Debug)]
pub struct SensitivityResult<S: Scalar> {
pub t: Vec<S>,
pub y: Vec<S>,
pub sensitivity: Vec<S>,
pub n_states: usize,
pub n_params: usize,
pub stats: SolverStats,
pub success: bool,
pub message: String,
}
impl<S: Scalar> SensitivityResult<S> {
pub fn len(&self) -> usize {
self.t.len()
}
pub fn is_empty(&self) -> bool {
self.t.is_empty()
}
pub fn y_at(&self, i: usize) -> &[S] {
let n = self.n_states;
let off = i * n;
&self.y[off..off + n]
}
pub fn sensitivity_at(&self, i: usize) -> &[S] {
let block = self.n_states * self.n_params;
let off = i * block;
&self.sensitivity[off..off + block]
}
pub fn sensitivity_for_param(&self, i: usize, k: usize) -> &[S] {
let n = self.n_states;
let block = n * self.n_params;
let off = i * block + k * n;
&self.sensitivity[off..off + n]
}
pub fn dyi_dpj(&self, i: usize, state: usize, param: usize) -> S {
let n = self.n_states;
let block = n * self.n_params;
self.sensitivity[i * block + param * n + state]
}
pub fn final_state(&self) -> &[S] {
if self.is_empty() {
&[]
} else {
self.y_at(self.len() - 1)
}
}
pub fn final_sensitivity(&self) -> &[S] {
if self.is_empty() {
&[]
} else {
self.sensitivity_at(self.len() - 1)
}
}
pub fn normalized_sensitivity_at(&self, i: usize, p_nominal: &[S]) -> Vec<S> {
assert_eq!(
p_nominal.len(),
self.n_params,
"p_nominal length {} does not match n_params {}",
p_nominal.len(),
self.n_params,
);
let n = self.n_states;
let np = self.n_params;
let y = self.y_at(i);
let s = self.sensitivity_at(i);
let mut y_max = S::ZERO;
for &v in y {
let av = v.abs();
if av > y_max {
y_max = av;
}
}
let zero_threshold = S::from_f64(1e-15) * (S::ONE + y_max);
let mut out = vec![S::ZERO; n * np];
for k in 0..np {
for j in 0..n {
let yj = y[j];
if yj.abs() <= zero_threshold {
out[k * n + j] = S::ZERO;
} else {
out[k * n + j] = (p_nominal[k] / yj) * s[k * n + j];
}
}
}
out
}
}
pub struct ClosureSystem<S: Scalar, F> {
rhs: F,
params: Vec<S>,
n_states: usize,
}
impl<S: Scalar, F> ClosureSystem<S, F>
where
F: Fn(S, &[S], &[S], &mut [S]),
{
pub fn new(rhs: F, params: Vec<S>, n_states: usize) -> Self {
Self {
rhs,
params,
n_states,
}
}
}
impl<S: Scalar, F> ParametricOdeSystem<S> for ClosureSystem<S, F>
where
F: Fn(S, &[S], &[S], &mut [S]),
{
fn n_states(&self) -> usize {
self.n_states
}
fn n_params(&self) -> usize {
self.params.len()
}
fn params(&self) -> &[S] {
&self.params
}
fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]) {
(self.rhs)(t, y, p, dydt)
}
}
pub fn solve_forward_sensitivity<Sol, S, Sys>(
system: &Sys,
t0: S,
tf: S,
y0: &[S],
options: &SolverOptions<S>,
) -> Result<SensitivityResult<S>, SolverError>
where
Sol: Solver<S>,
S: Scalar,
Sys: ParametricOdeSystem<S>,
{
let n = system.n_states();
let np = system.n_params();
assert_eq!(
y0.len(),
n,
"solve_forward_sensitivity: y0.len() = {} but n_states = {}",
y0.len(),
n,
);
assert_eq!(
system.params().len(),
np,
"solve_forward_sensitivity: params().len() = {} but n_params = {}",
system.params().len(),
np,
);
let aug = AugmentedSystem::new(system);
let z0 = aug.initial_augmented(y0);
let aug_dim = aug.augmented_dim();
let aug_result = Sol::solve(&aug, t0, tf, &z0, options)?;
if !aug_result.success {
return Ok(SensitivityResult {
t: Vec::new(),
y: Vec::new(),
sensitivity: Vec::new(),
n_states: n,
n_params: np,
stats: aug_result.stats,
success: false,
message: aug_result.message,
});
}
let n_times = aug_result.len();
let mut t = Vec::with_capacity(n_times);
let mut y = Vec::with_capacity(n_times * n);
let mut sensitivity = Vec::with_capacity(n_times * n * np);
for i in 0..n_times {
t.push(aug_result.t[i]);
let block_start = i * aug_dim;
let state_block = &aug_result.y[block_start..block_start + n];
let sens_block = &aug_result.y[block_start + n..block_start + aug_dim];
y.extend_from_slice(state_block);
sensitivity.extend_from_slice(sens_block);
}
Ok(SensitivityResult {
t,
y,
sensitivity,
n_states: n,
n_params: np,
stats: aug_result.stats,
success: true,
message: String::new(),
})
}
pub fn solve_forward_sensitivity_with<Sol, S, F>(
rhs: F,
y0: &[S],
params: &[S],
t0: S,
tf: S,
options: &SolverOptions<S>,
) -> Result<SensitivityResult<S>, SolverError>
where
Sol: Solver<S>,
S: Scalar,
F: Fn(S, &[S], &[S], &mut [S]),
{
let system = ClosureSystem {
rhs,
params: params.to_vec(),
n_states: y0.len(),
};
solve_forward_sensitivity::<Sol, S, _>(&system, t0, tf, y0, options)
}
#[cfg(test)]
mod tests {
use super::*;
struct ExpDecay {
k: f64,
}
impl ParametricOdeSystem<f64> for ExpDecay {
fn n_states(&self) -> usize {
1
}
fn n_params(&self) -> usize {
1
}
fn params(&self) -> &[f64] {
std::slice::from_ref(&self.k)
}
fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
dy[0] = -p[0] * y[0];
}
fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) {
jy[0] = -self.k;
}
fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) {
jp[0] = -y[0];
}
fn has_analytical_jacobian_y(&self) -> bool {
true
}
fn has_analytical_jacobian_p(&self) -> bool {
true
}
}
struct Lin2 {
p: [f64; 2],
}
impl ParametricOdeSystem<f64> for Lin2 {
fn n_states(&self) -> usize {
2
}
fn n_params(&self) -> usize {
2
}
fn params(&self) -> &[f64] {
&self.p
}
fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
dy[0] = -p[0] * y[0] + p[1] * y[1];
dy[1] = -p[1] * y[1];
}
fn jacobian_y(&self, _t: f64, _y: &[f64], jy: &mut [f64]) {
jy[0 * 2 + 0] = -self.p[0];
jy[0 * 2 + 1] = self.p[1];
jy[1 * 2 + 0] = 0.0;
jy[1 * 2 + 1] = -self.p[1];
}
fn jacobian_p(&self, _t: f64, y: &[f64], jp: &mut [f64]) {
jp[0 * 2 + 0] = -y[0];
jp[0 * 2 + 1] = 0.0;
jp[1 * 2 + 0] = y[1];
jp[1 * 2 + 1] = -y[1];
}
fn has_analytical_jacobian_y(&self) -> bool {
true
}
fn has_analytical_jacobian_p(&self) -> bool {
true
}
}
#[test]
fn augmented_dim_and_initial_state() {
let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
assert_eq!(aug.augmented_dim(), 2 * (1 + 2));
let z0 = aug.initial_augmented(&[1.0, 2.0]);
assert_eq!(z0.len(), 6);
assert_eq!(&z0[..2], &[1.0, 2.0]);
for s in &z0[2..] {
assert_eq!(*s, 0.0);
}
}
#[test]
fn augmented_rhs_matches_hand_derivation() {
let aug = AugmentedSystem::new(ExpDecay { k: 0.5 });
let z = vec![2.0, 1.0]; let mut dz = vec![0.0; 2];
aug.rhs(0.0, &z, &mut dz);
assert!((dz[0] + 1.0).abs() < 1e-12);
assert!((dz[1] + 2.5).abs() < 1e-12);
}
#[test]
fn augmented_rhs_matches_hand_derivation_two_param() {
let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
let z = vec![1.0, 2.0, 1.0, 0.0, 0.0, 1.0];
let mut dz = vec![0.0; 6];
aug.rhs(0.0, &z, &mut dz);
assert!((dz[0] - 0.0).abs() < 1e-12);
assert!((dz[1] + 1.0).abs() < 1e-12);
assert!((dz[2] + 2.0).abs() < 1e-12);
assert!((dz[3] - 0.0).abs() < 1e-12);
assert!((dz[4] - 2.5).abs() < 1e-12);
assert!((dz[5] + 2.5).abs() < 1e-12);
}
#[test]
fn augmented_jacobian_is_block_diagonal_with_jy() {
let aug = AugmentedSystem::new(Lin2 { p: [1.0, 0.5] });
let z = vec![1.0, 2.0, 1.0, 0.0, 0.0, 1.0];
let dim = aug.augmented_dim();
let mut jac = vec![0.0; dim * dim];
aug.jacobian(0.0, &z, &mut jac);
let jy = [-1.0_f64, 0.5, 0.0, -0.5];
for block in 0..3 {
let row0 = block * 2;
let col0 = block * 2;
for i in 0..2 {
for j in 0..2 {
let got = jac[(row0 + i) * dim + (col0 + j)];
let want = jy[i * 2 + j];
assert!(
(got - want).abs() < 1e-12,
"diag block {block}, ({i},{j}): got {got}, want {want}"
);
}
}
}
for br in 0..3 {
for bc in 0..3 {
if br == bc {
continue;
}
for i in 0..2 {
for j in 0..2 {
let r = br * 2 + i;
let c = bc * 2 + j;
assert_eq!(jac[r * dim + c], 0.0, "off-diag ({br},{bc}) ({i},{j})");
}
}
}
}
}
#[test]
fn fd_default_jacobian_y_matches_analytical() {
struct Lin2Fd {
p: [f64; 2],
}
impl ParametricOdeSystem<f64> for Lin2Fd {
fn n_states(&self) -> usize {
2
}
fn n_params(&self) -> usize {
2
}
fn params(&self) -> &[f64] {
&self.p
}
fn rhs_with_params(&self, _t: f64, y: &[f64], p: &[f64], dy: &mut [f64]) {
dy[0] = -p[0] * y[0] + p[1] * y[1];
dy[1] = -p[1] * y[1];
}
}
let analytical = Lin2 { p: [1.0, 0.5] };
let fd = Lin2Fd { p: [1.0, 0.5] };
let y = [1.0, 2.0];
let mut a_jy = [0.0; 4];
let mut f_jy = [0.0; 4];
analytical.jacobian_y(0.0, &y, &mut a_jy);
fd.jacobian_y(0.0, &y, &mut f_jy);
for k in 0..4 {
assert!(
(a_jy[k] - f_jy[k]).abs() < 1e-7,
"j_y[{k}] analytical={} fd={}",
a_jy[k],
f_jy[k]
);
}
let mut a_jp = [0.0; 4];
let mut f_jp = [0.0; 4];
analytical.jacobian_p(0.0, &y, &mut a_jp);
fd.jacobian_p(0.0, &y, &mut f_jp);
for k in 0..4 {
assert!(
(a_jp[k] - f_jp[k]).abs() < 1e-7,
"j_p[{k}] analytical={} fd={}",
a_jp[k],
f_jp[k]
);
}
}
#[test]
fn sensitivity_result_accessors() {
let res = SensitivityResult {
t: vec![0.0, 1.0],
y: vec![1.0, 2.0, 0.5, 1.5],
sensitivity: vec![
0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6,
],
n_states: 2,
n_params: 3,
stats: SolverStats::new(),
success: true,
message: String::new(),
};
assert_eq!(res.len(), 2);
assert_eq!(res.y_at(1), &[0.5, 1.5]);
assert_eq!(res.sensitivity_at(0).len(), 6);
assert_eq!(res.sensitivity_for_param(1, 2), &[1.5, 1.6]);
assert!((res.dyi_dpj(0, 1, 0) - 0.2).abs() < 1e-12);
assert!((res.dyi_dpj(1, 0, 2) - 1.5).abs() < 1e-12);
assert_eq!(res.final_state(), &[0.5, 1.5]);
assert_eq!(res.final_sensitivity().len(), 6);
}
#[test]
fn normalized_sensitivity_handles_zero_state_safely() {
let res = SensitivityResult {
t: vec![0.0],
y: vec![1.0, 0.0],
sensitivity: vec![0.5, 0.7],
n_states: 2,
n_params: 1,
stats: SolverStats::new(),
success: true,
message: String::new(),
};
let norm = res.normalized_sensitivity_at(0, &[2.0]);
assert!((norm[0] - 1.0).abs() < 1e-12);
assert_eq!(norm[1], 0.0);
}
use crate::DoPri5;
use crate::SolverOptions;
#[test]
fn solve_forward_sensitivity_with_dopri5_matches_analytical_decay() {
let result = solve_forward_sensitivity_with::<DoPri5, f64, _>(
|_t, y, p, dy| {
dy[0] = -p[0] * y[0];
},
&[1.0],
&[0.5],
0.0,
2.0,
&SolverOptions::default().rtol(1e-8).atol(1e-10),
)
.expect("solve failed");
assert!(result.success);
assert_eq!(result.n_states, 1);
assert_eq!(result.n_params, 1);
let last = result.len() - 1;
let dy_dk = result.dyi_dpj(last, 0, 0);
let t_last = result.t[last];
let analytical = -t_last * (-0.5 * t_last).exp();
assert!(
(dy_dk - analytical).abs() < 1e-5,
"computed {dy_dk}, analytical {analytical}, |err| = {}",
(dy_dk - analytical).abs()
);
}
#[test]
fn solve_trait_form_matches_closure_form() {
let trait_result = solve_forward_sensitivity::<DoPri5, _, _>(
&ExpDecay { k: 0.5 },
0.0,
2.0,
&[1.0],
&SolverOptions::default().rtol(1e-9).atol(1e-12),
)
.expect("trait solve failed");
let closure_result = solve_forward_sensitivity_with::<DoPri5, f64, _>(
|_t, y, p, dy| {
dy[0] = -p[0] * y[0];
},
&[1.0],
&[0.5],
0.0,
2.0,
&SolverOptions::default().rtol(1e-9).atol(1e-12),
)
.expect("closure solve failed");
let n_t = trait_result.len();
let n_c = closure_result.len();
let trait_final = trait_result.dyi_dpj(n_t - 1, 0, 0);
let closure_final = closure_result.dyi_dpj(n_c - 1, 0, 0);
assert!(
(trait_final - closure_final).abs() < 1e-5,
"trait {trait_final} vs closure {closure_final}, |err| = {}",
(trait_final - closure_final).abs()
);
}
#[test]
fn solve_propagates_solver_failure_as_unsuccessful_result() {
let opts = SolverOptions::default().rtol(1e-9).atol(1e-12).max_steps(1);
let outcome = solve_forward_sensitivity_with::<DoPri5, f64, _>(
|_t, y, p, dy| {
dy[0] = -p[0] * y[0];
},
&[1.0],
&[0.5],
0.0,
10.0,
&opts,
);
match outcome {
Ok(r) => {
assert!(!r.success, "expected unsuccessful result, got success");
assert!(r.t.is_empty());
assert!(r.y.is_empty());
assert!(r.sensitivity.is_empty());
}
Err(_) => {
}
}
}
}