use nereids_core::constants::{PIVOT_FLOOR, POISSON_EPSILON};
use crate::error::FittingError;
use crate::lm::{FitModel, FlatMatrix};
use crate::parameters::ParameterSet;
pub struct JointPoissonObjective<'a> {
pub model: &'a dyn FitModel,
pub o: &'a [f64],
pub s: &'a [f64],
pub c: f64,
pub active_mask: Option<&'a [bool]>,
}
impl<'a> JointPoissonObjective<'a> {
pub fn n_data(&self) -> usize {
self.o.len()
}
pub fn n_active(&self) -> usize {
crate::active_mask::active_count(self.active_mask, self.o.len())
}
#[inline]
fn bin_active(&self, i: usize) -> bool {
self.active_mask.is_none_or(|m| m[i])
}
fn validate_inputs(&self, t_len: usize) -> Result<(), FittingError> {
if self.s.len() != self.o.len() {
return Err(FittingError::LengthMismatch {
expected: self.o.len(),
actual: self.s.len(),
field: "sample_counts",
});
}
if let Some(m) = self.active_mask
&& m.len() != self.o.len()
{
return Err(FittingError::LengthMismatch {
expected: self.o.len(),
actual: m.len(),
field: "active_mask",
});
}
if !self.c.is_finite() || self.c <= 0.0 {
return Err(FittingError::InvalidConfig(format!(
"proton-charge ratio c = Q_s/Q_ob must be finite and > 0, got {}",
self.c
)));
}
if t_len != self.o.len() {
return Err(FittingError::LengthMismatch {
expected: self.o.len(),
actual: t_len,
field: "transmission",
});
}
validate_counts(self.o, "open_beam_counts")?;
validate_counts(self.s, "sample_counts")?;
Ok(())
}
#[inline]
pub fn profile_lambda(&self, t_i: f64, o_i: f64, s_i: f64) -> f64 {
let denom = 1.0 + self.c * t_i;
if denom <= POISSON_EPSILON {
0.0
} else {
self.c * (o_i + s_i) / denom
}
}
pub fn profile_lambda_per_bin(&self, t: &[f64]) -> Result<Vec<f64>, FittingError> {
self.validate_inputs(t.len())?;
Ok(t.iter()
.zip(self.o.iter())
.zip(self.s.iter())
.map(|((&ti, &oi), &si)| self.profile_lambda(ti, oi, si))
.collect())
}
pub fn deviance_from_transmission(&self, t: &[f64]) -> Result<f64, FittingError> {
self.validate_inputs(t.len())?;
let mut d = 0.0;
for (i, ((&t_i, &o_i), &s_i)) in t.iter().zip(self.o.iter()).zip(self.s.iter()).enumerate()
{
if !self.bin_active(i) {
continue;
}
d += binomial_deviance_term(s_i, o_i, t_i, self.c);
}
Ok(d)
}
pub fn deviance(&self, params: &[f64]) -> Result<f64, FittingError> {
let t = self.model.evaluate(params)?;
if t.len() != self.o.len() {
return Err(FittingError::LengthMismatch {
expected: self.o.len(),
actual: t.len(),
field: "transmission",
});
}
self.deviance_from_transmission(&t)
}
pub fn deviance_gradient_analytical(
&self,
params: &[f64],
free_param_indices: &[usize],
) -> Result<Option<Vec<f64>>, FittingError> {
let t = self.model.evaluate(params)?;
self.validate_inputs(t.len())?;
let jac = match self
.model
.analytical_jacobian(params, free_param_indices, &t)
{
Some(j) => j,
None => return Ok(None),
};
let n_free = free_param_indices.len();
let mut grad = vec![0.0f64; n_free];
for (i, (&t_i, (&o_i, &s_i))) in t.iter().zip(self.o.iter().zip(self.s.iter())).enumerate()
{
if !self.bin_active(i) {
continue;
}
let w = deviance_weight(s_i, o_i, t_i, self.c);
if w == 0.0 {
continue;
}
for (g, col) in grad.iter_mut().zip(0..n_free) {
let j = jac.get(i, col);
if j.is_finite() {
*g += w * j;
}
}
}
Ok(Some(grad))
}
pub fn fisher_information(
&self,
params: &[f64],
free_param_indices: &[usize],
) -> Result<Option<FlatMatrix>, FittingError> {
let t = self.model.evaluate(params)?;
self.validate_inputs(t.len())?;
let jac = match self
.model
.analytical_jacobian(params, free_param_indices, &t)
{
Some(j) => j,
None => return Ok(None),
};
let n_free = free_param_indices.len();
let mut info = FlatMatrix::zeros(n_free, n_free);
for (i, ((&t_i, &o_i), &s_i)) in t.iter().zip(self.o.iter()).zip(self.s.iter()).enumerate()
{
if !self.bin_active(i) {
continue;
}
let h = deviance_curvature(s_i, o_i, t_i, self.c);
if h == 0.0 {
continue;
}
for j in 0..n_free {
let jij = jac.get(i, j);
if !jij.is_finite() {
continue;
}
for k in 0..n_free {
let jik = jac.get(i, k);
if jik.is_finite() {
*info.get_mut(j, k) += h * jij * jik;
}
}
}
}
Ok(Some(info))
}
pub fn fisher_information_fd(
&self,
params: &mut ParameterSet,
fd_step: f64,
) -> Result<Option<FlatMatrix>, FittingError> {
let free_idx = params.free_indices();
let base_values = params.all_values();
let t_base = self.model.evaluate(&base_values)?;
self.validate_inputs(t_base.len())?;
let n_e = t_base.len();
let n_free = free_idx.len();
if n_free == 0 {
return Ok(Some(FlatMatrix::zeros(0, 0)));
}
let mut jac = FlatMatrix::zeros(n_e, n_free);
for (col, &idx) in free_idx.iter().enumerate() {
let original = params.params[idx].value;
let step = fd_step * (1.0 + original.abs());
params.params[idx].value = original + step;
params.params[idx].clamp();
let forward_step = params.params[idx].value - original;
let t_plus = if forward_step.abs() >= PIVOT_FLOOR {
Some(self.model.evaluate(¶ms.all_values())?)
} else {
None
};
params.params[idx].value = original - step;
params.params[idx].clamp();
let backward_step = original - params.params[idx].value;
let t_minus = if backward_step.abs() >= PIVOT_FLOOR {
Some(self.model.evaluate(¶ms.all_values())?)
} else {
None
};
params.params[idx].value = original;
let (t_a, t_b, denom) = match (t_plus, t_minus) {
(Some(tp), Some(tm)) => (tp, tm, forward_step + backward_step),
(Some(tp), None) => (tp, t_base.clone(), forward_step),
(None, Some(tm)) => (t_base.clone(), tm, backward_step),
(None, None) => continue,
};
if denom.abs() < PIVOT_FLOOR {
continue;
}
for i in 0..n_e {
let a = t_a[i];
let b = t_b[i];
if a.is_finite() && b.is_finite() {
*jac.get_mut(i, col) = (a - b) / denom;
}
}
}
let mut info = FlatMatrix::zeros(n_free, n_free);
for (i, ((&t_i, &o_i), &s_i)) in t_base
.iter()
.zip(self.o.iter())
.zip(self.s.iter())
.enumerate()
{
if !self.bin_active(i) {
continue;
}
let h = deviance_curvature(s_i, o_i, t_i, self.c);
if h == 0.0 {
continue;
}
for j in 0..n_free {
let jij = jac.get(i, j);
if !jij.is_finite() {
continue;
}
for k in 0..n_free {
let jik = jac.get(i, k);
if jik.is_finite() {
*info.get_mut(j, k) += h * jij * jik;
}
}
}
}
Ok(Some(info))
}
pub fn deviance_gradient_fd(
&self,
params: &mut ParameterSet,
fd_step: f64,
) -> Result<Vec<f64>, FittingError> {
let free_idx = params.free_indices();
let base_values = params.all_values();
let base_d = self.deviance(&base_values)?;
let mut grad = vec![0.0; free_idx.len()];
for (j, &idx) in free_idx.iter().enumerate() {
let original = params.params[idx].value;
let step = fd_step * (1.0 + original.abs());
params.params[idx].value = original + step;
params.params[idx].clamp();
let mut actual_step = params.params[idx].value - original;
if actual_step.abs() < PIVOT_FLOOR {
params.params[idx].value = original - step;
params.params[idx].clamp();
actual_step = params.params[idx].value - original;
if actual_step.abs() < PIVOT_FLOOR {
params.params[idx].value = original;
continue;
}
}
let perturbed_values = params.all_values();
let perturbed_d = match self.deviance(&perturbed_values) {
Ok(v) if v.is_finite() => v,
_ => {
params.params[idx].value = original;
continue;
}
};
params.params[idx].value = original;
grad[j] = (perturbed_d - base_d) / actual_step;
}
Ok(grad)
}
}
#[inline]
fn binomial_deviance_term(s: f64, o: f64, t: f64, c: f64) -> f64 {
debug_assert!(
s.is_finite() && s >= 0.0,
"binomial_deviance_term: S must be finite and >= 0, got {s}"
);
debug_assert!(
o.is_finite() && o >= 0.0,
"binomial_deviance_term: O must be finite and >= 0, got {o}"
);
debug_assert!(
c.is_finite() && c > 0.0,
"binomial_deviance_term: c must be finite and > 0, got {c}"
);
if !t.is_finite() {
return f64::NAN;
}
let t_safe = t.max(POISSON_EPSILON);
let n = s + o;
if n <= 0.0 {
return 0.0;
}
let ct = c * t_safe;
let one_plus_ct = 1.0 + ct;
let exp_s = ct / one_plus_ct * n; let exp_o = n / one_plus_ct;
let term_s = xlogy_ratio(s, exp_s);
let term_o = xlogy_ratio(o, exp_o);
2.0 * (term_s + term_o)
}
fn validate_counts(counts: &[f64], field: &'static str) -> Result<(), FittingError> {
for (i, &v) in counts.iter().enumerate() {
if !v.is_finite() || v < 0.0 {
return Err(FittingError::InvalidConfig(format!(
"{field}[{i}] must be finite and >= 0, got {v}"
)));
}
}
Ok(())
}
#[inline]
fn xlogy_ratio(x: f64, y: f64) -> f64 {
if x <= 0.0 {
0.0
} else {
let y_safe = y.max(POISSON_EPSILON);
x * (x / y_safe).ln()
}
}
#[inline]
fn deviance_weight(s: f64, o: f64, t: f64, c: f64) -> f64 {
if !t.is_finite() {
return 0.0;
}
let t_safe = t.max(POISSON_EPSILON);
let one_plus_ct = 1.0 + c * t_safe;
-2.0 * (s - o * c * t_safe) / (t_safe * one_plus_ct)
}
#[inline]
fn deviance_curvature(s: f64, o: f64, t: f64, c: f64) -> f64 {
if !t.is_finite() {
return 0.0;
}
let t_safe = t.max(POISSON_EPSILON);
let n = s + o;
let one_plus_ct = 1.0 + c * t_safe;
2.0 * n * c / (t_safe * one_plus_ct * one_plus_ct)
}
use crate::lm::{invert_matrix, solve_damped_system};
use crate::nelder_mead::{NelderMeadConfig, nelder_mead_minimize};
#[derive(Debug, Clone)]
pub struct JointPoissonFitConfig {
pub max_iter: usize,
pub lambda_init: f64,
pub lambda_up: f64,
pub lambda_down: f64,
pub armijo_c: f64,
pub backtrack: f64,
pub tol_d: f64,
pub tol_param: f64,
pub fd_step: f64,
pub enable_polish: bool,
pub polish: NelderMeadConfig,
pub compute_covariance: bool,
}
impl Default for JointPoissonFitConfig {
fn default() -> Self {
Self {
max_iter: 200,
lambda_init: 1e-3,
lambda_up: 10.0,
lambda_down: 0.1,
armijo_c: 1e-4,
backtrack: 0.5,
tol_d: 1e-8,
tol_param: 1e-8,
fd_step: 1e-6,
enable_polish: false,
polish: NelderMeadConfig {
xatol: 1e-9,
fatol: 1e-10,
max_iter: 5000,
initial_step_frac: 0.02,
initial_step_abs: 1e-4,
},
compute_covariance: true,
}
}
}
#[derive(Debug, Clone)]
pub struct JointPoissonResult {
pub deviance: f64,
pub deviance_per_dof: f64,
pub n_data: usize,
pub n_active: usize,
pub n_free: usize,
pub gn_iterations: usize,
pub polish_iterations: usize,
pub gn_converged: bool,
pub polish_converged: bool,
pub polish_improved: bool,
pub params: Vec<f64>,
pub covariance: Option<FlatMatrix>,
pub uncertainties: Option<Vec<f64>>,
}
pub fn joint_poisson_fit(
objective: &JointPoissonObjective<'_>,
params: &mut ParameterSet,
config: &JointPoissonFitConfig,
) -> Result<JointPoissonResult, FittingError> {
let n_data = objective.n_data();
if n_data == 0 {
return Err(FittingError::EmptyData);
}
if objective.s.len() != n_data {
return Err(FittingError::LengthMismatch {
expected: n_data,
actual: objective.s.len(),
field: "sample_counts",
});
}
if !objective.c.is_finite() || objective.c <= 0.0 {
return Err(FittingError::InvalidConfig(format!(
"proton-charge ratio c = Q_s/Q_ob must be finite and > 0, got {}",
objective.c
)));
}
validate_counts(objective.o, "open_beam_counts")?;
validate_counts(objective.s, "sample_counts")?;
if let Some(m) = objective.active_mask
&& m.len() != n_data
{
return Err(FittingError::LengthMismatch {
expected: n_data,
actual: m.len(),
field: "active_mask",
});
}
let n_free_initial = params.n_free();
let n_active_initial = objective.n_active();
if n_active_initial == 0 {
return Ok(JointPoissonResult {
deviance: f64::NAN,
deviance_per_dof: f64::NAN,
n_data,
n_active: 0,
n_free: n_free_initial,
gn_iterations: 0,
polish_iterations: 0,
gn_converged: false,
polish_converged: false,
polish_improved: false,
params: params.all_values(),
covariance: None,
uncertainties: None,
});
}
if n_active_initial < n_free_initial {
return Ok(JointPoissonResult {
deviance: f64::NAN,
deviance_per_dof: f64::NAN,
n_data,
n_active: n_active_initial,
n_free: n_free_initial,
gn_iterations: 0,
polish_iterations: 0,
gn_converged: false,
polish_converged: false,
polish_improved: false,
params: params.all_values(),
covariance: None,
uncertainties: None,
});
}
let stage1 = damped_fisher_stage(objective, params, config)?;
let best_d_stage1 = stage1.deviance;
let gn_iterations = stage1.iterations;
let gn_converged = stage1.converged;
let mut polish_iterations = 0usize;
let mut polish_converged = false;
let mut polish_improved = false;
let free_idx = params.free_indices();
if config.enable_polish && !free_idx.is_empty() && best_d_stage1.is_finite() {
let bounds: Vec<(f64, f64)> = free_idx
.iter()
.map(|&i| (params.params[i].lower, params.params[i].upper))
.collect();
let x0: Vec<f64> = free_idx.iter().map(|&i| params.params[i].value).collect();
let all_values_snapshot = params.all_values();
let obj_closure = |x: &[f64]| -> Result<f64, FittingError> {
let mut all = all_values_snapshot.clone();
for (j, &idx) in free_idx.iter().enumerate() {
all[idx] = x[j];
}
objective.deviance(&all)
};
let nm = nelder_mead_minimize(obj_closure, &x0, Some(&bounds), &config.polish)?;
polish_iterations = nm.iterations;
polish_converged = nm.self_converged;
if nm.fun < best_d_stage1 {
polish_improved = true;
for (j, &idx) in free_idx.iter().enumerate() {
params.params[idx].value = nm.x[j];
params.params[idx].clamp();
}
}
}
let final_values = params.all_values();
let final_deviance = objective.deviance(&final_values)?;
let n_free = params.n_free();
let n_active = objective.n_active();
let dof = n_active.saturating_sub(n_free);
let deviance_per_dof = if dof > 0 {
final_deviance / dof as f64
} else {
f64::NAN
};
let (covariance, uncertainties) = if config.compute_covariance {
let free_idx = params.free_indices();
let info_opt = match objective.fisher_information(&final_values, &free_idx)? {
Some(info) => Some(info),
None => objective.fisher_information_fd(params, config.fd_step)?,
};
match info_opt {
Some(info) => match invert_matrix(&info) {
Some(mut cov) => {
for v in cov.data.iter_mut() {
*v *= 2.0;
}
let u: Vec<f64> = (0..cov.nrows)
.map(|i| {
let v = cov.get(i, i);
if v > 0.0 { v.sqrt() } else { f64::NAN }
})
.collect();
(Some(cov), Some(u))
}
None => (None, None),
},
None => (None, None),
}
} else {
(None, None)
};
Ok(JointPoissonResult {
deviance: final_deviance,
deviance_per_dof,
n_data,
n_active,
n_free,
gn_iterations,
polish_iterations,
gn_converged,
polish_converged,
polish_improved,
params: final_values,
covariance,
uncertainties,
})
}
struct Stage1Output {
deviance: f64,
iterations: usize,
converged: bool,
}
fn damped_fisher_stage(
objective: &JointPoissonObjective<'_>,
params: &mut ParameterSet,
config: &JointPoissonFitConfig,
) -> Result<Stage1Output, FittingError> {
let mut lambda = config.lambda_init;
let mut iter = 0usize;
let mut converged = false;
let mut all_vals = params.all_values();
let mut d_current = objective.deviance(&all_vals)?;
while iter < config.max_iter {
iter += 1;
let free_idx = params.free_indices();
let n_free = free_idx.len();
if n_free == 0 {
converged = d_current.is_finite();
break;
}
let grad = match objective.deviance_gradient_analytical(&all_vals, &free_idx)? {
Some(g) => g,
None => objective.deviance_gradient_fd(params, config.fd_step)?,
};
let info = match objective.fisher_information(&all_vals, &free_idx)? {
Some(m) => m,
None => {
let mut ident = FlatMatrix::zeros(n_free, n_free);
for i in 0..n_free {
*ident.get_mut(i, i) = 1.0;
}
ident
}
};
let neg_grad: Vec<f64> = grad.iter().map(|&g| -g).collect();
let step = match solve_damped_system(&info, &neg_grad, lambda) {
Some(s) => s,
None => {
lambda *= config.lambda_up;
if lambda > 1e16 {
break;
}
continue;
}
};
let grad_dot_step = grad
.iter()
.zip(step.iter())
.map(|(&g, &s)| g * s)
.sum::<f64>();
let effective_step: Vec<f64> = if grad_dot_step >= 0.0 {
grad.iter().map(|&g| -g).collect()
} else {
step
};
let mut alpha = 1.0;
let mut accepted = false;
let d0 = d_current;
let mut trial_vals = all_vals.clone();
for _ in 0..50 {
for (j, &idx) in free_idx.iter().enumerate() {
trial_vals[idx] = all_vals[idx] + alpha * effective_step[j];
}
for &idx in free_idx.iter() {
let lo = params.params[idx].lower;
let hi = params.params[idx].upper;
if trial_vals[idx] < lo {
trial_vals[idx] = lo;
}
if trial_vals[idx] > hi {
trial_vals[idx] = hi;
}
}
let d_trial = match objective.deviance(&trial_vals) {
Ok(v) if v.is_finite() => v,
_ => f64::INFINITY,
};
let gdotp = grad
.iter()
.zip(effective_step.iter())
.map(|(&g, &s)| g * s)
.sum::<f64>();
if d_trial <= d0 + config.armijo_c * alpha * gdotp {
accepted = true;
break;
}
alpha *= config.backtrack;
if alpha < 1e-16 {
break;
}
}
if accepted {
for &idx in free_idx.iter() {
params.params[idx].value = trial_vals[idx];
params.params[idx].clamp();
}
let rel_change =
(d_current - objective.deviance(&trial_vals)?) / d_current.abs().max(1.0);
all_vals = params.all_values();
let new_d = objective.deviance(&all_vals)?;
let step_norm_sq = effective_step
.iter()
.map(|&s| (alpha * s).powi(2))
.sum::<f64>();
let step_norm = step_norm_sq.sqrt();
d_current = new_d;
lambda = (lambda * config.lambda_down).max(1e-16);
if rel_change.abs() < config.tol_d && step_norm < config.tol_param {
converged = true;
break;
}
} else {
lambda *= config.lambda_up;
if lambda > 1e16 {
break;
}
}
}
Ok(Stage1Output {
deviance: d_current,
iterations: iter,
converged,
})
}
#[cfg(test)]
mod tests {
use super::*;
use crate::parameters::FitParameter;
struct ConstModel {
n_e: usize,
}
impl FitModel for ConstModel {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(vec![params[0]; self.n_e])
}
fn analytical_jacobian(
&self,
_params: &[f64],
free_param_indices: &[usize],
y_current: &[f64],
) -> Option<FlatMatrix> {
let n_e = y_current.len();
let n_free = free_param_indices.len();
let mut jac = FlatMatrix::zeros(n_e, n_free);
for i in 0..n_e {
for (j, &pi) in free_param_indices.iter().enumerate() {
*jac.get_mut(i, j) = if pi == 0 { 1.0 } else { 0.0 };
}
}
Some(jac)
}
}
struct LinearModel<'a> {
e: &'a [f64],
}
impl<'a> FitModel for LinearModel<'a> {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(self
.e
.iter()
.map(|&ei| (params[0] - params[1] * ei).max(POISSON_EPSILON))
.collect())
}
fn analytical_jacobian(
&self,
_params: &[f64],
free_param_indices: &[usize],
y_current: &[f64],
) -> Option<FlatMatrix> {
let n_e = y_current.len();
let n_free = free_param_indices.len();
let mut jac = FlatMatrix::zeros(n_e, n_free);
for i in 0..n_e {
for (j, &pi) in free_param_indices.iter().enumerate() {
*jac.get_mut(i, j) = match pi {
0 => 1.0,
1 => -self.e[i],
_ => 0.0,
};
}
}
Some(jac)
}
}
#[test]
fn test_profile_lambda_closed_form_matches_bisection() {
let cases = [
(50.0_f64, 5.0_f64, 0.5_f64, 1.0_f64),
(1000.0, 900.0, 0.9, 5.98),
(10.0, 1.0, 0.1, 2.0),
(0.0, 5.0, 0.25, 1.5), (5.0, 0.0, 0.75, 3.0), ];
for (o, s, t, c) in cases {
let model = ConstModel { n_e: 1 };
let obj = JointPoissonObjective {
model: &model,
o: &[o],
s: &[s],
c,
active_mask: None,
};
let closed = obj.profile_lambda(t, o, s);
let score = |lam: f64| (o + s) / lam - (1.0 / c + t);
let (mut lo, mut hi) = (1e-10, 1e12);
assert!(score(lo) >= 0.0);
assert!(score(hi) <= 0.0);
for _ in 0..200 {
let mid = 0.5 * (lo + hi);
if score(mid) > 0.0 {
lo = mid;
} else {
hi = mid;
}
}
let bisect = 0.5 * (lo + hi);
let rel_err = ((closed - bisect) / bisect).abs();
assert!(
rel_err < 1e-9,
"profile λ̂ mismatch: closed={closed} bisect={bisect} rel_err={rel_err}"
);
}
}
#[test]
fn test_deviance_zero_at_exact_match() {
let t_val = 0.5;
let c = 2.0;
let n_bins = 5;
let o = vec![100.0; n_bins];
let s = vec![100.0; n_bins];
let t = vec![t_val; n_bins];
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c,
active_mask: None,
};
let d = obj.deviance_from_transmission(&t).unwrap();
assert!(d.abs() < 1e-8, "D should be ≈ 0 at exact match, got {d}");
let d_via_params = obj.deviance(&[t_val]).unwrap();
assert!(d_via_params.abs() < 1e-8);
}
#[test]
fn test_deviance_gradient_matches_fd() {
let e: Vec<f64> = (0..20).map(|i| 0.1 + 0.05 * i as f64).collect();
let theta_true = [0.95_f64, 0.1_f64];
let c = 3.0;
let lam = 500.0;
let model = LinearModel { e: &e };
let t_true = model.evaluate(&theta_true).unwrap();
let o: Vec<f64> = t_true.iter().map(|_| lam / c).collect();
let s: Vec<f64> = t_true.iter().map(|&ti| lam * ti).collect();
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c,
active_mask: None,
};
let theta_eval = [0.80_f64, 0.15_f64];
let free_idx = vec![0, 1];
let g_analytical = obj
.deviance_gradient_analytical(&theta_eval, &free_idx)
.unwrap()
.expect("LinearModel provides analytical jacobian");
let eps = 1e-6;
let mut g_fd = [0.0_f64; 2];
for j in 0..2 {
let mut tp = theta_eval;
let mut tm = theta_eval;
tp[j] += eps;
tm[j] -= eps;
let dp = obj.deviance(&tp).unwrap();
let dm = obj.deviance(&tm).unwrap();
g_fd[j] = (dp - dm) / (2.0 * eps);
}
for (a, f) in g_analytical.iter().zip(g_fd.iter()) {
let rel = ((a - f) / f.abs().max(1e-6)).abs();
assert!(
rel < 1e-4,
"analytical vs FD gradient disagree: analytical={a} fd={f} rel={rel}"
);
}
}
#[test]
fn test_deviance_per_dof_asymptote() {
let n_bins = 200;
let t_true = 0.35_f64;
let c = 2.0;
let lam = 50.0;
let n_reps = 30;
let mut d_per_dof_samples = Vec::with_capacity(n_reps);
let mut bias_samples = Vec::with_capacity(n_reps);
let mut rng = Xorshift(0xDEAD_BEEF_CAFE_BABE);
for _ in 0..n_reps {
let o: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam / c)).collect();
let s: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam * t_true)).collect();
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c,
active_mask: None,
};
let grid: Vec<f64> = (0..200).map(|i| 0.01 + 0.99 * (i as f64) / 199.0).collect();
let mut best = (grid[0], f64::INFINITY);
for &t_try in &grid {
let d_try = obj
.deviance_from_transmission(&vec![t_try; n_bins])
.unwrap();
if d_try < best.1 {
best = (t_try, d_try);
}
}
let dt = 0.01;
let (mut lo, mut hi) = ((best.0 - dt).max(POISSON_EPSILON), (best.0 + dt).min(0.999));
let grad_at = |t: f64| -> f64 {
let tvec = vec![t; n_bins];
let free_idx = [0_usize];
let g = obj
.deviance_gradient_analytical(&[t], &free_idx)
.unwrap()
.unwrap();
let _ = tvec; g[0]
};
let mut glo = grad_at(lo);
let mut ghi = grad_at(hi);
if glo * ghi < 0.0 {
for _ in 0..80 {
let mid = 0.5 * (lo + hi);
let gmid = grad_at(mid);
if gmid * glo < 0.0 {
hi = mid;
ghi = gmid;
} else {
lo = mid;
glo = gmid;
}
}
}
let t_hat = 0.5 * (lo + hi);
let d_hat = obj
.deviance_from_transmission(&vec![t_hat; n_bins])
.unwrap();
let dof = (n_bins - 1) as f64;
d_per_dof_samples.push(d_hat / dof);
bias_samples.push((t_hat - t_true) / t_true);
}
let mean_dpd: f64 = d_per_dof_samples.iter().sum::<f64>() / d_per_dof_samples.len() as f64;
let mean_bias: f64 = bias_samples.iter().sum::<f64>() / bias_samples.len() as f64;
assert!(
(0.85..=1.15).contains(&mean_dpd),
"D/(n-k) asymptote out of band: mean={mean_dpd}"
);
assert!(
mean_bias.abs() < 0.02,
"density bias > 2%: mean={mean_bias}"
);
}
#[test]
fn test_zero_counts_contribute_zero() {
let model = ConstModel { n_e: 3 };
let obj = JointPoissonObjective {
model: &model,
o: &[0.0, 10.0, 5.0],
s: &[0.0, 5.0, 2.0],
c: 1.5,
active_mask: None,
};
let d_full = obj.deviance_from_transmission(&[0.6, 0.6, 0.6]).unwrap();
let obj_reduced = JointPoissonObjective {
model: &model, o: &[10.0, 5.0],
s: &[5.0, 2.0],
c: 1.5,
active_mask: None,
};
let d_reduced = obj_reduced.deviance_from_transmission(&[0.6, 0.6]).unwrap();
assert!((d_full - d_reduced).abs() < 1e-12);
}
#[test]
fn test_fd_gradient_matches_analytical() {
let e: Vec<f64> = (0..15).map(|i| 0.2 + 0.1 * i as f64).collect();
let theta = [0.9_f64, 0.05_f64];
let c = 1.5;
let lam = 300.0;
let model = LinearModel { e: &e };
let t_true = model.evaluate(&theta).unwrap();
let o: Vec<f64> = t_true.iter().map(|_| lam / c).collect();
let s: Vec<f64> = t_true.iter().map(|&ti| lam * ti).collect();
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c,
active_mask: None,
};
let mut ps = ParameterSet::new(vec![
FitParameter::non_negative("theta_0", 0.85),
FitParameter::non_negative("theta_1", 0.06),
]);
let g_fd = obj.deviance_gradient_fd(&mut ps, 1e-6).unwrap();
let g_analytical = obj
.deviance_gradient_analytical(&ps.all_values(), &ps.free_indices())
.unwrap()
.unwrap();
for (f, a) in g_fd.iter().zip(g_analytical.iter()) {
let rel = ((f - a) / a.abs().max(1e-6)).abs();
assert!(rel < 5e-3, "fd={f} analytical={a} rel={rel}");
}
}
#[test]
fn test_fisher_matrix_symmetry_psd() {
let e: Vec<f64> = (0..10).map(|i| 0.3 + 0.1 * i as f64).collect();
let theta = [0.9_f64, 0.05_f64];
let c = 2.0;
let lam = 400.0;
let model = LinearModel { e: &e };
let t_true = model.evaluate(&theta).unwrap();
let o: Vec<f64> = t_true.iter().map(|_| lam / c).collect();
let s: Vec<f64> = t_true.iter().map(|&ti| lam * ti).collect();
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c,
active_mask: None,
};
let info = obj
.fisher_information(&theta, &[0, 1])
.unwrap()
.expect("LinearModel provides analytical jacobian");
let i01 = info.get(0, 1);
let i10 = info.get(1, 0);
assert!((i01 - i10).abs() < 1e-10);
assert!(info.get(0, 0) > 0.0);
assert!(info.get(1, 1) > 0.0);
let det = info.get(0, 0) * info.get(1, 1) - i01 * i10;
assert!(det > 0.0, "Fisher matrix determinant = {det}");
}
struct BackgroundedTransmission<'a> {
inner: &'a dyn FitModel,
energies: &'a [f64],
n_idx: usize,
a_idx: usize,
b_a_idx: usize,
b_b_idx: usize,
b_c_idx: usize,
n_params: usize,
}
impl<'a> FitModel for BackgroundedTransmission<'a> {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let t_inner = self.inner.evaluate(&[params[self.n_idx]])?;
let a_n = params[self.a_idx];
let b_a = params[self.b_a_idx];
let b_b = params[self.b_b_idx];
let b_c = params[self.b_c_idx];
Ok(t_inner
.iter()
.zip(self.energies.iter())
.map(|(&t, &e)| {
let inv_sqrt_e = if e > 0.0 { 1.0 / e.sqrt() } else { 0.0 };
let sqrt_e = if e > 0.0 { e.sqrt() } else { 0.0 };
a_n * t + b_a + b_b * inv_sqrt_e + b_c * sqrt_e
})
.collect())
}
}
struct ExpDecayModel<'a> {
sigma: &'a [f64],
}
impl<'a> FitModel for ExpDecayModel<'a> {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let n = params[0];
Ok(self
.sigma
.iter()
.map(|&s| (-n * s).exp().max(POISSON_EPSILON))
.collect())
}
fn analytical_jacobian(
&self,
_params: &[f64],
free_param_indices: &[usize],
y_current: &[f64],
) -> Option<FlatMatrix> {
let n_e = y_current.len();
let n_free = free_param_indices.len();
let mut jac = FlatMatrix::zeros(n_e, n_free);
for (i, &y_i) in y_current.iter().enumerate() {
for (j, &pi) in free_param_indices.iter().enumerate() {
*jac.get_mut(i, j) = if pi == 0 { -self.sigma[i] * y_i } else { 0.0 };
}
}
Some(jac)
}
}
struct Xorshift(u64);
impl Xorshift {
fn next_u64(&mut self) -> u64 {
let mut x = self.0;
x ^= x << 13;
x ^= x >> 7;
x ^= x << 17;
self.0 = x;
x
}
fn uniform(&mut self) -> f64 {
(self.next_u64() as f64) / (u64::MAX as f64)
}
fn poisson(&mut self, lambda: f64) -> f64 {
if lambda <= 0.0 {
return 0.0;
}
if lambda > 30.0 {
let u1 = self.uniform().max(1e-12);
let u2 = self.uniform();
let z = (-2.0 * u1.ln()).sqrt() * (2.0 * std::f64::consts::PI * u2).cos();
return (lambda + z * lambda.sqrt()).round().max(0.0);
}
let l = (-lambda).exp();
let mut k: f64 = 0.0;
let mut p: f64 = 1.0;
loop {
k += 1.0;
let u = self.uniform();
p *= u;
if p <= l {
return k - 1.0;
}
if k > 1000.0 {
return k - 1.0;
}
}
}
}
#[test]
fn test_joint_poisson_fit_matched_model_single_param() {
let n_bins = 200;
let sigma = vec![1.0_f64; n_bins];
let model = ExpDecayModel { sigma: &sigma };
let n_true = 0.3_f64;
let c = 5.98;
let lam = 3000.0; let t_true = model.evaluate(&[n_true]).unwrap();
let mut rng = Xorshift(0x1234_5678_9ABC_DEF0);
let o: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam / c)).collect();
let s: Vec<f64> = (0..n_bins).map(|i| rng.poisson(lam * t_true[i])).collect();
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("n", 0.1)]);
let cfg = JointPoissonFitConfig {
enable_polish: true,
..Default::default()
};
let result = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
let n_fit = result.params[0];
let rel_bias = (n_fit - n_true) / n_true;
assert!(
rel_bias.abs() < 0.01,
"density bias {rel_bias} exceeds 1% (n_fit={n_fit} n_true={n_true})"
);
assert!(
(0.85..=1.15).contains(&result.deviance_per_dof),
"D/(n-k) out of band: {}",
result.deviance_per_dof
);
}
#[test]
fn test_joint_poisson_fit_polish_does_not_worsen_deviance() {
let n_bins = 150;
let energies: Vec<f64> = (0..n_bins).map(|i| 1.0 + 0.5 * i as f64).collect();
let sigma: Vec<f64> = energies.iter().map(|&e| 1.0 / e).collect();
let inner = ExpDecayModel { sigma: &sigma };
let n_true = 0.3_f64;
let a_n_true = 0.9_f64;
let t_inner_true = inner.evaluate(&[n_true]).unwrap();
let t_true: Vec<f64> = t_inner_true.iter().map(|&t| a_n_true * t).collect();
let c = 5.98_f64;
let lam = 5000.0_f64;
let mut rng = Xorshift(0xF00D_FACE_DEAD_BEEF);
let o: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam / c)).collect();
let s: Vec<f64> = (0..n_bins).map(|i| rng.poisson(lam * t_true[i])).collect();
let bg_model = BackgroundedTransmission {
inner: &inner,
energies: &energies,
n_idx: 0,
a_idx: 1,
b_a_idx: 2,
b_b_idx: 3,
b_c_idx: 4,
n_params: 5,
};
let _ = bg_model.n_params;
let obj = JointPoissonObjective {
model: &bg_model,
o: &o,
s: &s,
c,
active_mask: None,
};
let mk_params = || {
ParameterSet::new(vec![
FitParameter::non_negative("n", 0.25),
FitParameter::non_negative("A_n", 1.0),
FitParameter {
name: "B_A".into(),
value: 0.0,
lower: -0.05,
upper: 0.05,
fixed: false,
},
FitParameter {
name: "B_B".into(),
value: 0.0,
lower: -0.05,
upper: 0.05,
fixed: false,
},
FitParameter {
name: "B_C".into(),
value: 0.0,
lower: -0.05,
upper: 0.05,
fixed: false,
},
])
};
let mut params_no_polish = mk_params();
let cfg_no_polish = JointPoissonFitConfig {
enable_polish: false,
..Default::default()
};
let r_no_polish = joint_poisson_fit(&obj, &mut params_no_polish, &cfg_no_polish).unwrap();
let mut params_polish = mk_params();
let cfg_polish = JointPoissonFitConfig {
enable_polish: true,
..Default::default()
};
let r_polish = joint_poisson_fit(&obj, &mut params_polish, &cfg_polish).unwrap();
assert!(
r_polish.deviance <= r_no_polish.deviance + 1e-6,
"polish worsened D: D_polish={} D_no_polish={}",
r_polish.deviance,
r_no_polish.deviance
);
if r_polish.polish_improved {
assert!(
r_polish.deviance < r_no_polish.deviance,
"polish_improved=true but D_polish={} >= D_no_polish={}",
r_polish.deviance,
r_no_polish.deviance
);
}
let n_fit = r_polish.params[0];
assert!(n_fit.is_finite() && n_fit > 0.0);
assert!(
n_fit > 0.1 && n_fit < 0.8,
"density grossly off: n_fit={n_fit} (truth={n_true})"
);
}
#[test]
fn test_joint_poisson_fit_exposes_separate_converged_flags() {
let n_bins = 50;
let sigma = vec![0.5_f64; n_bins];
let model = ExpDecayModel { sigma: &sigma };
let n_true = 0.2;
let c = 2.0;
let lam = 500.0;
let t_true = model.evaluate(&[n_true]).unwrap();
let mut rng = Xorshift(0xABAD_CAFE_BABE_F00D);
let o: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam / c)).collect();
let s: Vec<f64> = (0..n_bins).map(|i| rng.poisson(lam * t_true[i])).collect();
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("n", 0.1)]);
let cfg = JointPoissonFitConfig {
enable_polish: true,
..Default::default()
};
let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
assert!(r.gn_converged || r.polish_converged);
assert!(r.n_data == n_bins);
assert!(r.n_free == 1);
assert!(r.deviance > 0.0);
assert!(r.deviance_per_dof.is_finite());
assert!(r.uncertainties.is_some());
let u = r.uncertainties.as_ref().unwrap();
assert_eq!(u.len(), 1);
assert!(u[0].is_finite() && u[0] > 0.0);
}
#[test]
fn test_uncertainty_matches_analytical_fisher_inverse() {
let n_bins = 200;
let t_true = 0.5_f64;
let c = 2.0_f64;
let lam = 100.0_f64;
let o: Vec<f64> = vec![lam / c; n_bins];
let s: Vec<f64> = vec![lam * t_true; n_bins];
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", t_true)]);
let cfg = JointPoissonFitConfig {
enable_polish: false,
..Default::default()
};
let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
let sigma_reported = r.uncertainties.as_ref().expect("σ available")[0];
let sigma_analytical = (t_true * (1.0 + c * t_true) / (n_bins as f64 * lam)).sqrt();
let rel_err = (sigma_reported - sigma_analytical).abs() / sigma_analytical;
assert!(
rel_err < 0.05,
"reported σ = {sigma_reported} vs analytical I^{{-1}}^(1/2) = \
{sigma_analytical} (rel_err = {rel_err}); pre-fix code reported \
σ_analytical / √2 ≈ {} which would give rel_err ≈ 0.293",
sigma_analytical / 2.0_f64.sqrt(),
);
}
#[test]
fn test_jp_active_mask_subset_equivalence() {
let o_full = [10.0, 20.0, 5.0, 15.0, 25.0];
let s_full = [4.0, 8.0, 1.0, 6.0, 12.0];
let t_full = [0.4, 0.5, 0.7, 0.6, 0.45];
let mask = [false, true, false, true, false];
let c = 1.5;
let model_full = ConstModel { n_e: 5 };
let obj_full = JointPoissonObjective {
model: &model_full,
o: &o_full,
s: &s_full,
c,
active_mask: Some(&mask),
};
let d_masked = obj_full.deviance_from_transmission(&t_full).unwrap();
let o_sub = [o_full[1], o_full[3]];
let s_sub = [s_full[1], s_full[3]];
let t_sub = [t_full[1], t_full[3]];
let model_sub = ConstModel { n_e: 2 };
let obj_sub = JointPoissonObjective {
model: &model_sub,
o: &o_sub,
s: &s_sub,
c,
active_mask: None,
};
let d_subset = obj_sub.deviance_from_transmission(&t_sub).unwrap();
assert!(
(d_masked - d_subset).abs() < 1e-12,
"masked deviance {d_masked} != subset deviance {d_subset}"
);
assert_eq!(obj_full.n_active(), 2);
assert_eq!(obj_full.n_data(), 5);
}
#[test]
fn test_jp_active_mask_gradient_subset_equivalence() {
let e_full: Vec<f64> = (0..6).map(|i| 0.1 + 0.1 * i as f64).collect();
let theta_true = [0.95_f64, 0.05_f64];
let c = 2.0;
let lam = 100.0;
let model_full = LinearModel { e: &e_full };
let t_full = model_full.evaluate(&theta_true).unwrap();
let o_full: Vec<f64> = vec![lam / c; e_full.len()];
let s_full: Vec<f64> = t_full.iter().map(|&ti| lam * ti).collect();
let mask = vec![false, false, true, true, true, false];
let obj_full = JointPoissonObjective {
model: &model_full,
o: &o_full,
s: &s_full,
c,
active_mask: Some(&mask),
};
let params_full = ParameterSet::new(vec![
FitParameter::non_negative("a", theta_true[0]),
FitParameter::non_negative("b", theta_true[1]),
]);
let free_idx = params_full.free_indices();
let theta_eval = [0.9_f64, 0.07_f64];
let grad_masked = obj_full
.deviance_gradient_analytical(&theta_eval, &free_idx)
.unwrap()
.expect("analytical gradient");
let e_sub = e_full[2..5].to_vec();
let o_sub = o_full[2..5].to_vec();
let s_sub = s_full[2..5].to_vec();
let model_sub = LinearModel { e: &e_sub };
let obj_sub = JointPoissonObjective {
model: &model_sub,
o: &o_sub,
s: &s_sub,
c,
active_mask: None,
};
let grad_sub = obj_sub
.deviance_gradient_analytical(&theta_eval, &free_idx)
.unwrap()
.expect("analytical gradient");
for (i, (&gm, &gs)) in grad_masked.iter().zip(grad_sub.iter()).enumerate() {
assert!(
(gm - gs).abs() < 1e-9,
"grad component {i}: masked={gm} subset={gs}"
);
}
}
#[test]
fn test_joint_poisson_rejects_zero_active_mask() {
let n_bins = 10;
let o: Vec<f64> = vec![50.0; n_bins];
let s: Vec<f64> = vec![25.0; n_bins];
let mask = vec![false; n_bins]; let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: Some(&mask),
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
let cfg = JointPoissonFitConfig::default();
let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
assert!(
!r.gn_converged && !r.polish_converged,
"underdetermined fit must report non-converged"
);
assert!(
r.deviance.is_nan(),
"underdetermined deviance must be NaN, got {}",
r.deviance
);
assert!(
r.deviance_per_dof.is_nan(),
"underdetermined deviance-per-dof must be NaN, got {}",
r.deviance_per_dof
);
assert_eq!(r.n_data, n_bins);
assert_eq!(r.n_active, 0);
assert_eq!(r.n_free, 1);
assert!(r.covariance.is_none());
assert!(r.uncertainties.is_none());
}
#[test]
fn test_joint_poisson_rejects_zero_active_with_no_free_params() {
let n_bins = 5;
let o: Vec<f64> = vec![10.0; n_bins];
let s: Vec<f64> = vec![5.0; n_bins];
let mask = vec![false; n_bins];
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: Some(&mask),
};
let mut params = ParameterSet::new(vec![FitParameter::fixed("T", 0.5)]);
let r = joint_poisson_fit(&obj, &mut params, &JointPoissonFitConfig::default()).unwrap();
assert!(!r.gn_converged);
assert!(!r.polish_converged);
assert!(r.deviance.is_nan());
assert!(r.deviance_per_dof.is_nan());
assert_eq!(r.n_active, 0);
assert_eq!(r.n_free, 0);
}
#[test]
fn test_joint_poisson_rejects_active_mask_length_mismatch() {
let n_bins = 5;
let o: Vec<f64> = vec![10.0; n_bins];
let s: Vec<f64> = vec![5.0; n_bins];
let mask_wrong = vec![true, true, true]; let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: Some(&mask_wrong),
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
let cfg = JointPoissonFitConfig::default();
let err = joint_poisson_fit(&obj, &mut params, &cfg).unwrap_err();
assert!(
matches!(
err,
FittingError::LengthMismatch {
field: "active_mask",
..
}
),
"expected LengthMismatch on active_mask; got {err:?}"
);
}
#[test]
fn test_joint_poisson_rejects_o_s_length_mismatch() {
let n_bins = 5;
let o: Vec<f64> = vec![10.0; n_bins];
let s: Vec<f64> = vec![5.0; n_bins - 1];
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
let err =
joint_poisson_fit(&obj, &mut params, &JointPoissonFitConfig::default()).unwrap_err();
assert!(
matches!(
err,
FittingError::LengthMismatch {
field: "sample_counts",
..
}
),
"expected LengthMismatch on sample_counts; got {err:?}"
);
}
#[test]
fn test_joint_poisson_rejects_non_positive_c() {
let n_bins = 5;
let o: Vec<f64> = vec![10.0; n_bins];
let s: Vec<f64> = vec![5.0; n_bins];
let model = ConstModel { n_e: n_bins };
let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
let obj_zero = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 0.0,
active_mask: None,
};
let err = joint_poisson_fit(&obj_zero, &mut params, &JointPoissonFitConfig::default())
.unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(_)),
"expected InvalidConfig on c=0; got {err:?}"
);
let mut params2 = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
let obj_neg = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: -1.5,
active_mask: None,
};
let err = joint_poisson_fit(&obj_neg, &mut params2, &JointPoissonFitConfig::default())
.unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(_)),
"expected InvalidConfig on c<0; got {err:?}"
);
let mut params3 = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
let obj_nan = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: f64::NAN,
active_mask: None,
};
let err = joint_poisson_fit(&obj_nan, &mut params3, &JointPoissonFitConfig::default())
.unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(_)),
"expected InvalidConfig on c=NaN; got {err:?}"
);
}
#[test]
fn test_binomial_deviance_term_nan_t_returns_nan() {
let d_nan_t = binomial_deviance_term(50.0, 10.0, f64::NAN, 2.0);
assert!(
d_nan_t.is_nan(),
"non-finite T must produce NaN deviance, not a finite shim; got {d_nan_t}"
);
let d_inf_t = binomial_deviance_term(50.0, 10.0, f64::INFINITY, 2.0);
assert!(
d_inf_t.is_nan(),
"+inf T must produce NaN deviance; got {d_inf_t}"
);
let d_neg_inf_t = binomial_deviance_term(50.0, 10.0, f64::NEG_INFINITY, 2.0);
assert!(
d_neg_inf_t.is_nan(),
"-inf T must produce NaN deviance; got {d_neg_inf_t}"
);
}
#[test]
fn test_deviance_weight_nan_t_returns_zero() {
let w = deviance_weight(50.0, 10.0, f64::NAN, 2.0);
assert_eq!(w, 0.0, "non-finite T must give zero weight; got {w}");
}
#[test]
fn test_deviance_curvature_nan_t_returns_zero() {
let h = deviance_curvature(50.0, 10.0, f64::NAN, 2.0);
assert_eq!(h, 0.0, "non-finite T must give zero curvature; got {h}");
}
#[test]
fn test_joint_poisson_fit_rejects_nan_transmission() {
struct NanAtSmallTheta;
impl FitModel for NanAtSmallTheta {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let t = if params[0] < 0.1 { f64::NAN } else { 0.5 };
Ok(vec![t; 4])
}
fn analytical_jacobian(
&self,
_params: &[f64],
free_param_indices: &[usize],
y_current: &[f64],
) -> Option<FlatMatrix> {
let n_e = y_current.len();
let n_free = free_param_indices.len();
let mut jac = FlatMatrix::zeros(n_e, n_free);
for i in 0..n_e {
for (j, &pi) in free_param_indices.iter().enumerate() {
*jac.get_mut(i, j) = if pi == 0 { 1.0 } else { 0.0 };
}
}
Some(jac)
}
}
let model = NanAtSmallTheta;
let n = 4;
let o = vec![10.0; n];
let s = vec![5.0; n];
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.05)]);
let cfg = JointPoissonFitConfig::default();
let result = joint_poisson_fit(&obj, &mut params, &cfg);
match result {
Ok(r) => {
if r.params[0] < 0.1 {
assert!(
r.deviance.is_nan() && !r.gn_converged,
"stayed in NaN region but reported finite deviance: {r:?}"
);
}
}
Err(_) => {
}
}
}
#[test]
fn test_joint_poisson_all_fixed_nan_transmission_does_not_converge() {
struct NanModel {
n_e: usize,
}
impl FitModel for NanModel {
fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(vec![f64::NAN; self.n_e])
}
}
let n_bins = 5;
let o = vec![10.0; n_bins];
let s = vec![5.0; n_bins];
let model = NanModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::fixed("T", 0.5)]);
let cfg = JointPoissonFitConfig::default();
let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
assert!(
r.deviance.is_nan(),
"expected NaN deviance from all-fixed NaN model; got {}",
r.deviance
);
assert!(
r.deviance_per_dof.is_nan(),
"expected NaN deviance_per_dof; got {}",
r.deviance_per_dof
);
assert!(
!r.gn_converged,
"all-fixed NaN deviance must not be reported as GN-converged",
);
assert_eq!(r.n_free, 0);
assert_eq!(r.n_active, n_bins);
assert_eq!(
r.gn_iterations, 1,
"all-fixed branch should report exactly one iteration",
);
}
#[test]
fn test_joint_poisson_all_fixed_nan_transmission_with_polish_does_not_panic() {
struct NanModel {
n_e: usize,
}
impl FitModel for NanModel {
fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(vec![f64::NAN; self.n_e])
}
}
let n_bins = 5;
let o = vec![10.0; n_bins];
let s = vec![5.0; n_bins];
let model = NanModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::fixed("T", 0.5)]);
let cfg = JointPoissonFitConfig {
enable_polish: true,
..JointPoissonFitConfig::default()
};
let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
assert!(
r.deviance.is_nan(),
"expected NaN deviance from all-fixed NaN model; got {}",
r.deviance
);
assert!(
!r.gn_converged,
"all-fixed NaN deviance must not be reported as GN-converged",
);
assert!(
!r.polish_converged,
"polish stage must report not-converged when skipped on all-fixed params",
);
assert!(
!r.polish_improved,
"polish stage cannot have improved the deviance when it was skipped",
);
assert_eq!(
r.polish_iterations, 0,
"polish stage must report zero iterations when skipped",
);
assert_eq!(r.n_free, 0);
assert_eq!(r.n_active, n_bins);
assert_eq!(
r.gn_iterations, 1,
"all-fixed branch should report exactly one iteration",
);
}
#[test]
fn test_joint_poisson_polish_does_not_report_converged_when_stage1_nan() {
struct NanModel {
n_e: usize,
}
impl FitModel for NanModel {
fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(vec![f64::NAN; self.n_e])
}
}
let n_bins = 5;
let o = vec![10.0; n_bins];
let s = vec![5.0; n_bins];
let model = NanModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
let cfg = JointPoissonFitConfig {
enable_polish: true,
..JointPoissonFitConfig::default()
};
let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
assert!(
r.deviance.is_nan(),
"expected NaN deviance from NaN model; got {}",
r.deviance
);
assert!(!r.gn_converged, "stage 1 cannot converge on NaN deviance",);
assert!(
!r.polish_converged,
"stage 2 must not report converged when stage 1 ended non-finite",
);
assert!(
!r.polish_improved,
"polish cannot have improved a NaN starting deviance",
);
assert_eq!(
r.polish_iterations, 0,
"polish must not run when stage 1 is non-finite",
);
assert_eq!(r.n_free, 1);
assert_eq!(r.n_active, n_bins);
}
#[test]
fn test_fisher_information_fd_skips_nan_probe() {
struct NanFdProbe;
impl FitModel for NanFdProbe {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let t = if (params[0] - 0.6).abs() > 1e-3 {
f64::NAN
} else {
params[0]
};
Ok(vec![t; 3])
}
}
let model = NanFdProbe;
let n = 3;
let o = vec![10.0; n];
let s = vec![5.0; n];
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.6)]);
let info = obj
.fisher_information_fd(&mut params, 1e-2)
.expect("fisher_information_fd should not return Err on a finite base")
.expect("fisher_information_fd should return Some(matrix)");
for v in info.data.iter() {
assert!(
v.is_finite(),
"fisher_information_fd produced non-finite entry: {v}"
);
}
}
#[test]
fn test_deviance_from_transmission_rejects_non_finite_counts() {
let n_bins = 4;
let mut o = vec![10.0; n_bins];
o[2] = f64::NAN;
let s = vec![5.0; n_bins];
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let t = vec![0.5; n_bins];
let err = obj.deviance_from_transmission(&t).unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(ref msg) if msg.contains("open_beam_counts")),
"expected InvalidConfig naming open_beam_counts; got {err:?}"
);
let mut s_inf = vec![5.0; n_bins];
s_inf[0] = f64::INFINITY;
let obj_inf = JointPoissonObjective {
model: &model,
o: &vec![10.0; n_bins],
s: &s_inf,
c: 1.0,
active_mask: None,
};
let err = obj_inf.deviance_from_transmission(&t).unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(ref msg) if msg.contains("sample_counts")),
"expected InvalidConfig naming sample_counts; got {err:?}"
);
}
#[test]
fn test_deviance_from_transmission_rejects_negative_counts() {
let n_bins = 3;
let mut o = vec![10.0; n_bins];
o[1] = -2.0;
let s = vec![5.0; n_bins];
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let t = vec![0.5; n_bins];
let err = obj.deviance_from_transmission(&t).unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(ref msg) if msg.contains("open_beam_counts")),
"expected InvalidConfig naming open_beam_counts; got {err:?}"
);
}
#[test]
fn test_other_public_methods_reject_non_finite_counts() {
let n_bins = 4;
let mut s = vec![5.0; n_bins];
s[3] = f64::NAN;
let o = vec![10.0; n_bins];
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let t = vec![0.5; n_bins];
let err = obj.profile_lambda_per_bin(&t).unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(_)),
"profile_lambda_per_bin: expected InvalidConfig; got {err:?}"
);
let params = vec![0.5];
let free_idx = vec![0];
let err = obj
.deviance_gradient_analytical(¶ms, &free_idx)
.unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(_)),
"deviance_gradient_analytical: expected InvalidConfig; got {err:?}"
);
let err = obj.fisher_information(¶ms, &free_idx).unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(_)),
"fisher_information: expected InvalidConfig; got {err:?}"
);
let mut ps = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
let err = obj.fisher_information_fd(&mut ps, 1e-2).unwrap_err();
assert!(
matches!(err, FittingError::InvalidConfig(_)),
"fisher_information_fd: expected InvalidConfig; got {err:?}"
);
}
#[test]
fn test_validate_inputs_reports_transmission_length_mismatch_correctly() {
let n_bins = 5;
let o = vec![10.0; n_bins];
let s = vec![5.0; n_bins];
let model = ConstModel { n_e: n_bins };
let obj = JointPoissonObjective {
model: &model,
o: &o,
s: &s,
c: 1.0,
active_mask: None,
};
let t_short = vec![0.5; n_bins - 2];
let err = obj.deviance_from_transmission(&t_short).unwrap_err();
match err {
FittingError::LengthMismatch {
expected,
actual,
field,
} => {
assert_eq!(field, "transmission", "field must name `transmission`");
assert_eq!(expected, n_bins, "expected must be o.len()");
assert_eq!(actual, n_bins - 2, "actual must be t.len()");
}
other => panic!("expected LengthMismatch on transmission; got {other:?}"),
}
}
}