use numra_core::uncertainty::Uncertain;
use numra_core::Scalar;
use crate::error::SolverError;
use crate::problem::OdeSystem;
use crate::sensitivity::solve_forward_sensitivity_with;
use crate::solver::{Solver, SolverOptions, SolverResult, SolverStats};
#[derive(Clone, Debug)]
pub enum UncertaintyMode {
Trajectory,
MonteCarlo {
n_samples: usize,
},
}
#[derive(Clone, Debug)]
pub struct UncertainParam<S: Scalar> {
pub name: String,
pub nominal: S,
pub std: S,
}
impl<S: Scalar> UncertainParam<S> {
pub fn new(name: impl Into<String>, nominal: S, std: S) -> Self {
Self {
name: name.into(),
nominal,
std,
}
}
pub fn from_uncertain(name: impl Into<String>, u: Uncertain<S>) -> Self {
Self {
name: name.into(),
nominal: u.mean,
std: u.std(),
}
}
pub fn variance(&self) -> S {
self.std * self.std
}
}
#[derive(Clone, Debug)]
pub struct UncertainSolverResult<S: Scalar> {
pub result: SolverResult<S>,
pub sigma: Vec<S>,
pub sensitivities: Option<Vec<Vec<S>>>,
pub params: Vec<UncertainParam<S>>,
}
impl<S: Scalar> UncertainSolverResult<S> {
pub fn sigma_at(&self, i: usize, j: usize) -> S {
self.sigma[i * self.result.dim + j]
}
pub fn uncertain_at(&self, i: usize, j: usize) -> Uncertain<S> {
let mean = self.result.y_at(i)[j];
let std = self.sigma_at(i, j);
Uncertain::from_std(mean, std)
}
pub fn sensitivity_at(&self, i: usize, j: usize, k: usize) -> Option<S> {
self.sensitivities.as_ref().map(|sens| {
let n_params = self.params.len();
sens[i][j * n_params + k]
})
}
pub fn len(&self) -> usize {
self.result.len()
}
pub fn is_empty(&self) -> bool {
self.result.is_empty()
}
}
pub fn solve_trajectory<Sol, S, F>(
model: F,
y0: &[S],
t0: S,
tf: S,
params: &[UncertainParam<S>],
options: &SolverOptions<S>,
) -> Result<UncertainSolverResult<S>, SolverError>
where
S: Scalar,
Sol: Solver<S>,
F: Fn(S, &[S], &mut [S], &[S]),
{
let n_states = y0.len();
let n_params = params.len();
let nominal_params: Vec<S> = params.iter().map(|p| p.nominal).collect();
let variances: Vec<S> = params.iter().map(|p| p.variance()).collect();
let rhs = move |t: S, y: &[S], p: &[S], dydt: &mut [S]| {
model(t, y, dydt, p);
};
let sens =
solve_forward_sensitivity_with::<Sol, S, _>(rhs, y0, &nominal_params, t0, tf, options)?;
if !sens.success {
return Err(SolverError::Other(sens.message));
}
let n_times = sens.len();
let mut sens_out = Vec::with_capacity(n_times);
let mut sigma_out = Vec::with_capacity(n_times * n_states);
for i in 0..n_times {
let block = sens.sensitivity_at(i);
let mut row_major = vec![S::ZERO; n_states * n_params];
for j in 0..n_states {
for k in 0..n_params {
row_major[j * n_params + k] = block[k * n_states + j];
}
}
sens_out.push(row_major);
for j in 0..n_states {
let mut var_j = S::ZERO;
for k in 0..n_params {
let dydp = block[k * n_states + j];
var_j = var_j + dydp * dydp * variances[k];
}
sigma_out.push(var_j.sqrt());
}
}
let nominal_result = SolverResult {
t: sens.t,
y: sens.y,
dim: n_states,
stats: sens.stats,
success: true,
message: String::new(),
events: Vec::new(),
terminated_by_event: false,
dense_output: None,
};
Ok(UncertainSolverResult {
result: nominal_result,
sigma: sigma_out,
sensitivities: Some(sens_out),
params: params.to_vec(),
})
}
pub fn solve_monte_carlo<Sol, S, F>(
model: F,
y0: &[S],
t0: S,
tf: S,
params: &[UncertainParam<S>],
n_samples: usize,
options: &SolverOptions<S>,
seed: u64,
) -> Result<UncertainSolverResult<S>, SolverError>
where
S: Scalar,
Sol: Solver<S>,
F: Fn(S, &[S], &mut [S], &[S]) + Send + Sync,
{
let n_states = y0.len();
let n_params = params.len();
let nominal_params: Vec<S> = params.iter().map(|p| p.nominal).collect();
let nominal_sys = ParameterizedWrapper {
model: &model,
params: nominal_params.clone(),
n_dim: n_states,
};
let nominal_result = Sol::solve(&nominal_sys, t0, tf, y0, options)?;
if !nominal_result.success {
return Err(SolverError::Other(nominal_result.message));
}
let n_times = nominal_result.len();
let mut sum_final = vec![S::ZERO; n_states];
let mut sum_sq_final = vec![S::ZERO; n_states];
let mut n_success: usize = 0;
let mut rng_state = seed;
for _ in 0..n_samples {
let mut p_sample = Vec::with_capacity(n_params);
for param in params {
let z = box_muller_sample(&mut rng_state);
let p_val = param.nominal + param.std * S::from_f64(z);
p_sample.push(p_val);
}
let sample_sys = ParameterizedWrapper {
model: &model,
params: p_sample,
n_dim: n_states,
};
match Sol::solve(&sample_sys, t0, tf, y0, options) {
Ok(result) if result.success => {
if let Some(y_final) = result.y_final() {
n_success += 1;
for j in 0..n_states {
sum_final[j] = sum_final[j] + y_final[j];
sum_sq_final[j] = sum_sq_final[j] + y_final[j] * y_final[j];
}
}
}
_ => {
}
}
}
if n_success < 2 {
return Err(SolverError::Other(
"Monte Carlo: fewer than 2 samples succeeded".to_string(),
));
}
let n_s = S::from_usize(n_success);
let mut sigma_final = Vec::with_capacity(n_states);
for j in 0..n_states {
let mean = sum_final[j] / n_s;
let var = (sum_sq_final[j] / n_s - mean * mean) * n_s / (n_s - S::ONE);
let std = if var > S::ZERO { var.sqrt() } else { S::ZERO };
sigma_final.push(std);
}
let mut sigma = Vec::with_capacity(n_times * n_states);
for i in 0..n_times {
let frac = if n_times > 1 {
S::from_usize(i) / S::from_usize(n_times - 1)
} else {
S::ONE
};
for j in 0..n_states {
sigma.push(sigma_final[j] * frac);
}
}
let mc_result = SolverResult {
t: nominal_result.t.clone(),
y: nominal_result.y.clone(),
dim: n_states,
stats: SolverStats::new(),
success: true,
message: format!("{}/{} samples succeeded", n_success, n_samples),
events: Vec::new(),
terminated_by_event: false,
dense_output: None,
};
Ok(UncertainSolverResult {
result: mc_result,
sigma,
sensitivities: None,
params: params.to_vec(),
})
}
pub fn solve_with_uncertainty<Sol, S, F>(
model: F,
y0: &[S],
t0: S,
tf: S,
params: &[UncertainParam<S>],
mode: &UncertaintyMode,
options: &SolverOptions<S>,
seed: Option<u64>,
) -> Result<UncertainSolverResult<S>, SolverError>
where
S: Scalar,
Sol: Solver<S>,
F: Fn(S, &[S], &mut [S], &[S]) + Send + Sync,
{
match mode {
UncertaintyMode::Trajectory => {
solve_trajectory::<Sol, S, F>(model, y0, t0, tf, params, options)
}
UncertaintyMode::MonteCarlo { n_samples } => solve_monte_carlo::<Sol, S, F>(
model,
y0,
t0,
tf,
params,
*n_samples,
options,
seed.unwrap_or(42),
),
}
}
struct ParameterizedWrapper<'a, S: Scalar, F> {
model: &'a F,
params: Vec<S>,
n_dim: usize,
}
impl<S: Scalar, F> OdeSystem<S> for ParameterizedWrapper<'_, S, F>
where
F: Fn(S, &[S], &mut [S], &[S]) + Send + Sync,
{
fn dim(&self) -> usize {
self.n_dim
}
fn rhs(&self, t: S, y: &[S], dydt: &mut [S]) {
(self.model)(t, y, dydt, &self.params);
}
}
fn splitmix64(state: &mut u64) -> u64 {
*state = state.wrapping_add(0x9e3779b97f4a7c15);
let mut z = *state;
z = (z ^ (z >> 30)).wrapping_mul(0xbf58476d1ce4e5b9);
z = (z ^ (z >> 27)).wrapping_mul(0x94d049bb133111eb);
z ^ (z >> 31)
}
fn box_muller_sample(state: &mut u64) -> f64 {
loop {
let u1 = ((splitmix64(state) >> 11) as f64) / ((1u64 << 53) as f64);
let u2 = ((splitmix64(state) >> 11) as f64) / ((1u64 << 53) as f64);
if u1 > 1e-15 {
let r = (-2.0 * u1.ln()).sqrt();
return r * (2.0 * core::f64::consts::PI * u2).cos();
}
}
}
#[cfg(test)]
mod tests {
use super::*;
use crate::{DoPri5, Radau5};
#[test]
fn test_trajectory_exponential_decay() {
let params = vec![UncertainParam::new("k", 0.5, 0.05)];
let result = solve_with_uncertainty::<DoPri5, f64, _>(
|_t, y, dydt, p| {
dydt[0] = -p[0] * y[0];
},
&[1.0],
0.0,
5.0,
¶ms,
&UncertaintyMode::Trajectory,
&SolverOptions::default().rtol(1e-8).atol(1e-10),
None,
)
.expect("solve_with_uncertainty failed");
assert!(result.result.success);
assert!(!result.is_empty());
let n = result.len();
let t_final = result.result.t[n - 1];
let k = 0.5;
let sigma_k = 0.05;
let exact_sigma = t_final * (-k * t_final).exp() * sigma_k;
let computed_sigma = result.sigma_at(n - 1, 0);
assert!(
(computed_sigma - exact_sigma).abs() < 0.001,
"sigma: computed={}, exact={}, err={}",
computed_sigma,
exact_sigma,
(computed_sigma - exact_sigma).abs()
);
assert!(result.sensitivities.is_some());
let dydp = result.sensitivity_at(n - 1, 0, 0).unwrap();
let exact_dydp = -t_final * (-k * t_final).exp();
assert!(
(dydp - exact_dydp).abs() < 0.001,
"dy/dk: computed={}, exact={}, err={}",
dydp,
exact_dydp,
(dydp - exact_dydp).abs()
);
}
#[test]
fn test_trajectory_two_params() {
let params = vec![
UncertainParam::new("a", 1.0, 0.1),
UncertainParam::new("b", 2.0, 0.2),
];
let result = solve_trajectory::<DoPri5, f64, _>(
|_t, y, dydt, p| {
dydt[0] = -p[0] * y[0] + p[1];
},
&[1.0],
0.0,
3.0,
¶ms,
&SolverOptions::default().rtol(1e-8).atol(1e-10),
)
.expect("solve_trajectory failed");
assert!(result.result.success);
let n = result.len();
let y_final = result.result.y_at(n - 1)[0];
assert!(
(y_final - 2.0).abs() < 0.1,
"y(3) = {}, expected near 2.0",
y_final
);
let sigma_final = result.sigma_at(n - 1, 0);
assert!(
sigma_final > 0.0,
"sigma should be positive: {}",
sigma_final
);
}
#[test]
fn test_trajectory_vs_monte_carlo() {
let params = vec![UncertainParam::new("k", 0.5, 0.05)];
let traj_result = solve_with_uncertainty::<DoPri5, f64, _>(
|_t, y, dydt, p| {
dydt[0] = -p[0] * y[0];
},
&[1.0],
0.0,
2.0,
¶ms,
&UncertaintyMode::Trajectory,
&SolverOptions::default().rtol(1e-8).atol(1e-10),
None,
)
.expect("trajectory failed");
let mc_result = solve_with_uncertainty::<DoPri5, f64, _>(
|_t, y, dydt, p| {
dydt[0] = -p[0] * y[0];
},
&[1.0],
0.0,
2.0,
¶ms,
&UncertaintyMode::MonteCarlo { n_samples: 5000 },
&SolverOptions::default().rtol(1e-8).atol(1e-10),
Some(12345),
)
.expect("monte carlo failed");
let n_traj = traj_result.len();
let n_mc = mc_result.len();
let sigma_traj = traj_result.sigma_at(n_traj - 1, 0);
let sigma_mc = mc_result.sigma_at(n_mc - 1, 0);
let rel_diff = (sigma_traj - sigma_mc).abs() / sigma_traj;
assert!(
rel_diff < 0.15,
"Trajectory sigma={}, MC sigma={}, rel_diff={}",
sigma_traj,
sigma_mc,
rel_diff
);
}
#[test]
fn test_composability_stiff_solver() {
let params = vec![UncertainParam::new("k", 50.0, 5.0)];
let result = solve_trajectory::<Radau5, f64, _>(
|_t, y, dydt, p| {
dydt[0] = -p[0] * y[0];
},
&[1.0],
0.0,
0.2,
¶ms,
&SolverOptions::default().rtol(1e-4).atol(1e-7).h0(1e-4),
)
.expect("stiff solve failed");
assert!(result.result.success);
assert!(!result.is_empty());
let n = result.len();
let sigma = result.sigma_at(n - 1, 0);
assert!(sigma >= 0.0, "sigma should be non-negative: {}", sigma);
}
#[test]
fn test_trajectory_lotka_volterra() {
let params = vec![
UncertainParam::new("alpha", 1.0, 0.1),
UncertainParam::new("beta", 0.1, 0.01),
UncertainParam::new("delta", 0.075, 0.005),
UncertainParam::new("gamma", 1.5, 0.1),
];
let result = solve_trajectory::<DoPri5, f64, _>(
|_t, y, dydt, p| {
let x = y[0];
let yy = y[1];
dydt[0] = p[0] * x - p[1] * x * yy;
dydt[1] = p[2] * x * yy - p[3] * yy;
},
&[10.0, 5.0],
0.0,
5.0,
¶ms,
&SolverOptions::default().rtol(1e-8).atol(1e-10),
)
.expect("Lotka-Volterra solve failed");
assert!(result.result.success);
assert_eq!(result.result.dim, 2);
let n = result.len();
let sigma_prey = result.sigma_at(n - 1, 0);
let sigma_pred = result.sigma_at(n - 1, 1);
assert!(sigma_prey > 0.0, "prey sigma should be positive");
assert!(sigma_pred > 0.0, "predator sigma should be positive");
let sens = result.sensitivities.as_ref().unwrap();
assert_eq!(sens[n - 1].len(), 2 * 4);
}
#[test]
fn test_uncertain_param_from_uncertain() {
let u = Uncertain::from_std(5.0, 0.5);
let p = UncertainParam::from_uncertain("x", u);
assert!((p.nominal - 5.0).abs() < 1e-10);
assert!((p.std - 0.5).abs() < 1e-10);
assert!((p.variance() - 0.25).abs() < 1e-10);
}
#[test]
fn test_zero_uncertainty() {
let params = vec![UncertainParam::new("k", 0.5, 0.0)];
let result = solve_trajectory::<DoPri5, f64, _>(
|_t, y, dydt, p| {
dydt[0] = -p[0] * y[0];
},
&[1.0],
0.0,
2.0,
¶ms,
&SolverOptions::default(),
)
.expect("solve failed");
for i in 0..result.len() {
let sigma = result.sigma_at(i, 0);
assert!(
sigma.abs() < 1e-10,
"sigma should be ~0 at t={}: got {}",
result.result.t[i],
sigma
);
}
}
#[test]
fn test_uncertain_at() {
let params = vec![UncertainParam::new("k", 0.5, 0.05)];
let result = solve_trajectory::<DoPri5, f64, _>(
|_t, y, dydt, p| {
dydt[0] = -p[0] * y[0];
},
&[1.0],
0.0,
1.0,
¶ms,
&SolverOptions::default().rtol(1e-8),
)
.expect("solve failed");
let n = result.len();
let u = result.uncertain_at(n - 1, 0);
assert!(u.mean > 0.0);
assert!(u.variance > 0.0);
assert!((u.std() - result.sigma_at(n - 1, 0)).abs() < 1e-14);
}
}