nereids_fitting/joint_poisson.rs
1//! Joint-Poisson counts-path objective with profiled flux.
2//!
3//! This module implements the **joint-Poisson conditional binomial deviance**
4//! derived in `.research/spatial-regularization/evidence/32-counts-path-governing-equations-v2.md`
5//! (equations §4.1, §5.7, §6.2b) and validated experimentally in memo 35
6//! §P1. It supersedes the fixed-flux Poisson NLL (`poisson.rs`) for the
7//! counts-path fitter.
8//!
9//! ## Model
10//!
11//! Under the λ-at-sample convention with proton-charge ratio `c = Q_s / Q_ob`:
12//!
13//! - `O_i ~ Poisson(λ_i / c)` (open-beam counts)
14//! - `S_i ~ Poisson(λ_i · T_i)` (sample counts)
15//!
16//! Profiling out `λ_i` bin-by-bin gives the closed-form MLE
17//!
18//! ```text
19//! λ̂_i = c · (O_i + S_i) / (1 + c · T_i)
20//! ```
21//!
22//! The profile-conditional log-likelihood is equivalent (up to constants) to
23//! a Binomial `S_i | N_i = O_i + S_i ~ Binomial(N_i, p_i)` with
24//!
25//! ```text
26//! p_i = c · T_i / (1 + c · T_i)
27//! ```
28//!
29//! The conditional deviance is
30//!
31//! ```text
32//! D(θ) = 2 · Σ_i [ S_i · ln(S_i / (N_i · p_i))
33//! + O_i · ln(O_i / (N_i · (1 − p_i))) ]
34//! ```
35//!
36//! with the `x · ln(x / 0) → 0` convention when `x = 0`.
37//!
38//! Under the correct model, `D / (n − k)` → 1 as n → ∞ — this replaces the
39//! fixed-flux Pearson χ²/dof reported from the old Poisson path (which
40//! scaled with `c` at constant density fidelity; see memo 35 headline).
41
42use nereids_core::constants::{PIVOT_FLOOR, POISSON_EPSILON};
43
44use crate::error::FittingError;
45use crate::lm::{FitModel, FlatMatrix};
46use crate::parameters::ParameterSet;
47
48/// Joint-Poisson objective.
49///
50/// Wraps a transmission `FitModel` (which produces `T_i = model.evaluate(θ)`)
51/// together with the observed open-beam counts `O_i`, sample counts `S_i`,
52/// and proton-charge ratio `c = Q_s / Q_ob`.
53///
54/// The caller is responsible for ensuring `o`, `s`, and `model.evaluate()`
55/// output all have the same length.
56pub struct JointPoissonObjective<'a> {
57 /// Transmission model: `evaluate(θ) → T(E)`.
58 pub model: &'a dyn FitModel,
59 /// Open-beam counts per bin.
60 pub o: &'a [f64],
61 /// Sample counts per bin.
62 pub s: &'a [f64],
63 /// Proton-charge ratio `c = Q_s / Q_ob`. Must be strictly positive.
64 pub c: f64,
65 /// Optional per-bin active mask (SAMMY EMIN/EMAX-equivalent
66 /// fit-energy-range restriction). When `Some(m)`, only bins where
67 /// `m[i]` is `true` contribute to the deviance / gradient / Fisher
68 /// information; the model is still evaluated on the full grid so
69 /// resolution broadening at the boundaries is correct. When
70 /// `None`, all bins are active (default behaviour).
71 ///
72 /// Length must equal `o.len()`; the GUI / pipeline dispatch builds
73 /// it from the configured `[E_min, E_max]` against the energy grid.
74 pub active_mask: Option<&'a [bool]>,
75}
76
77impl<'a> JointPoissonObjective<'a> {
78 /// Number of data bins.
79 pub fn n_data(&self) -> usize {
80 self.o.len()
81 }
82
83 /// Number of *active* data bins — `n_data` when no mask is set,
84 /// or the count of `true` entries in `active_mask` otherwise.
85 /// This is the count that should drive deviance-per-dof reporting.
86 pub fn n_active(&self) -> usize {
87 crate::active_mask::active_count(self.active_mask, self.o.len())
88 }
89
90 /// Predicate: is bin `i` active? Returns `true` when no mask is
91 /// set (full-grid default).
92 #[inline]
93 fn bin_active(&self, i: usize) -> bool {
94 self.active_mask.is_none_or(|m| m[i])
95 }
96
97 /// Runtime guard for the public methods that bypass `joint_poisson_fit`'s
98 /// up-front validation (callers may invoke `deviance_from_transmission`,
99 /// `deviance_gradient_analytical`, `fisher_information[_fd]`, etc.
100 /// directly for diagnostics). Mirrors the entry-point checks in
101 /// `joint_poisson_fit`: `o.len() == s.len()`, `c` finite and > 0, optional
102 /// `active_mask` length agrees, all `o[i]` / `s[i]` finite and >= 0, and
103 /// the caller-supplied transmission length agrees with `o.len()`. The
104 /// `debug_assert!`s in the per-bin helpers are no-ops in release builds —
105 /// without this guard a length mismatch in `s` would silently truncate
106 /// via `.zip()` and a non-positive / NaN `c` would produce finite
107 /// garbage.
108 ///
109 /// **Error orientation.** `FittingError::LengthMismatch` displays as
110 /// `"{field} length ({actual}) must match expected length ({expected})"`.
111 /// The objective's own invariants (`s.len()` vs `o.len()`, `mask.len()`
112 /// vs `o.len()`) are checked first, with `expected = o.len()` so the
113 /// message accurately names the offending field. The caller-supplied
114 /// `t` length is then checked against `o.len()` with `field =
115 /// "transmission"` — pre-fix this branch reported `field =
116 /// "open_beam_counts"` with `expected = t_len`, which read as "the
117 /// open-beam array is wrong" when the actual fault was the caller's
118 /// transmission slice.
119 fn validate_inputs(&self, t_len: usize) -> Result<(), FittingError> {
120 // Internal invariants of the objective itself — these must hold
121 // regardless of what the caller passes for `t`.
122 if self.s.len() != self.o.len() {
123 return Err(FittingError::LengthMismatch {
124 expected: self.o.len(),
125 actual: self.s.len(),
126 field: "sample_counts",
127 });
128 }
129 if let Some(m) = self.active_mask
130 && m.len() != self.o.len()
131 {
132 return Err(FittingError::LengthMismatch {
133 expected: self.o.len(),
134 actual: m.len(),
135 field: "active_mask",
136 });
137 }
138 if !self.c.is_finite() || self.c <= 0.0 {
139 return Err(FittingError::InvalidConfig(format!(
140 "proton-charge ratio c = Q_s/Q_ob must be finite and > 0, got {}",
141 self.c
142 )));
143 }
144 // Caller-supplied length: the transmission slice must match the
145 // objective's bin count.
146 if t_len != self.o.len() {
147 return Err(FittingError::LengthMismatch {
148 expected: self.o.len(),
149 actual: t_len,
150 field: "transmission",
151 });
152 }
153 // Per-element count validation. The entry point `joint_poisson_fit`
154 // also calls `validate_counts` up-front so the user gets the error
155 // before any LM work, but every public method that bypasses the
156 // entry point (`deviance_from_transmission`, `fisher_information`,
157 // `profile_lambda_per_bin`, …) must still reject non-finite /
158 // negative counts — the inner `binomial_deviance_term` /
159 // `xlogy_ratio` would otherwise propagate NaN past the zero-clamp
160 // (`NaN <= 0.0` is `false`) or silently swallow a negative as zero.
161 validate_counts(self.o, "open_beam_counts")?;
162 validate_counts(self.s, "sample_counts")?;
163 Ok(())
164 }
165
166 /// Closed-form profile MLE for the per-bin flux: `λ̂ = c·(O+S) / (1+c·T)`.
167 ///
168 /// Guards: when `1 + c·T ≤ ε`, returns 0 to avoid division blow-up.
169 #[inline]
170 pub fn profile_lambda(&self, t_i: f64, o_i: f64, s_i: f64) -> f64 {
171 let denom = 1.0 + self.c * t_i;
172 if denom <= POISSON_EPSILON {
173 0.0
174 } else {
175 self.c * (o_i + s_i) / denom
176 }
177 }
178
179 /// Vector form of [`profile_lambda`](Self::profile_lambda).
180 ///
181 /// Validates `t.len() == o.len() == s.len()` and `c > 0`; returns
182 /// `FittingError::LengthMismatch` / `InvalidConfig` rather than the
183 /// previous `.zip()` truncate-and-pretend behaviour (which would
184 /// silently shrink the output to `min(t.len(), o.len(), s.len())`).
185 pub fn profile_lambda_per_bin(&self, t: &[f64]) -> Result<Vec<f64>, FittingError> {
186 self.validate_inputs(t.len())?;
187 Ok(t.iter()
188 .zip(self.o.iter())
189 .zip(self.s.iter())
190 .map(|((&ti, &oi), &si)| self.profile_lambda(ti, oi, si))
191 .collect())
192 }
193
194 /// Conditional binomial deviance at the given transmission vector.
195 ///
196 /// D = 2 · Σ [ S·ln(S/(Np)) + O·ln(O/(N(1−p))) ] with
197 /// `p = cT/(1+cT)`, `N = O+S`, and `x·ln(x/0) → 0`.
198 ///
199 /// Near invalid or numerically tiny transmission values, the per-bin
200 /// evaluation (`binomial_deviance_term`) uses `t.max(POISSON_EPSILON)`
201 /// to clamp T away from zero before entering the logarithms and the
202 /// `1/(1+cT)` factor. This avoids singular logs and division-by-zero
203 /// but is a piecewise clamp, not a smooth quadratic extrapolation —
204 /// D(T) is C⁰ at the clamp boundary, not C¹. In practice this is
205 /// adequate because the optimizer's transmission values come from a
206 /// `FitModel` that keeps T bounded well above `POISSON_EPSILON` for
207 /// physically plausible density / nuisance parameter values.
208 pub fn deviance_from_transmission(&self, t: &[f64]) -> Result<f64, FittingError> {
209 self.validate_inputs(t.len())?;
210 let mut d = 0.0;
211 for (i, ((&t_i, &o_i), &s_i)) in t.iter().zip(self.o.iter()).zip(self.s.iter()).enumerate()
212 {
213 if !self.bin_active(i) {
214 continue;
215 }
216 d += binomial_deviance_term(s_i, o_i, t_i, self.c);
217 }
218 Ok(d)
219 }
220
221 /// Evaluate the deviance at parameter vector θ by calling the model.
222 pub fn deviance(&self, params: &[f64]) -> Result<f64, FittingError> {
223 let t = self.model.evaluate(params)?;
224 if t.len() != self.o.len() {
225 return Err(FittingError::LengthMismatch {
226 expected: self.o.len(),
227 actual: t.len(),
228 field: "transmission",
229 });
230 }
231 self.deviance_from_transmission(&t)
232 }
233
234 /// Analytical gradient of the deviance w.r.t. the free parameters.
235 ///
236 /// Returns `None` if the transmission model does not provide an analytical
237 /// Jacobian — callers should fall back to `deviance_gradient_fd`.
238 ///
239 /// Gradient derivation: with `p_i = cT_i/(1+cT_i)` and N_i = O_i+S_i,
240 ///
241 /// d D / d T_i = −2 · (S_i − O_i·c·T_i) / (T_i · (1 + c·T_i))
242 ///
243 /// then chain-rule with the transmission Jacobian J_{i,j} = ∂T_i / ∂θ_{f(j)}
244 /// where f(j) is the j-th free parameter index.
245 pub fn deviance_gradient_analytical(
246 &self,
247 params: &[f64],
248 free_param_indices: &[usize],
249 ) -> Result<Option<Vec<f64>>, FittingError> {
250 let t = self.model.evaluate(params)?;
251 self.validate_inputs(t.len())?;
252 let jac = match self
253 .model
254 .analytical_jacobian(params, free_param_indices, &t)
255 {
256 Some(j) => j,
257 None => return Ok(None),
258 };
259 let n_free = free_param_indices.len();
260 let mut grad = vec![0.0f64; n_free];
261 for (i, (&t_i, (&o_i, &s_i))) in t.iter().zip(self.o.iter().zip(self.s.iter())).enumerate()
262 {
263 if !self.bin_active(i) {
264 continue;
265 }
266 let w = deviance_weight(s_i, o_i, t_i, self.c);
267 // `deviance_weight` returns 0 for non-finite `t_i`, so a NaN
268 // transmission row already contributes nothing — except that
269 // `0.0 * NaN = NaN`. If the upstream Jacobian column has a
270 // NaN cell (common for FD-built Jacobians where the model
271 // returns NaN at some probe point), the bare `0.0 * jac.get(...)`
272 // would poison `grad[col]`. Skip the row entirely when the
273 // weight is zero, and skip any individual Jacobian cell that
274 // is not finite.
275 if w == 0.0 {
276 continue;
277 }
278 for (g, col) in grad.iter_mut().zip(0..n_free) {
279 let j = jac.get(i, col);
280 if j.is_finite() {
281 *g += w * j;
282 }
283 }
284 }
285 Ok(Some(grad))
286 }
287
288 /// Fisher information for free parameters (Gauss-Newton curvature of D).
289 ///
290 /// Uses the expected-info form
291 ///
292 /// h_i ≡ ∂² D / ∂ T_i² ≈ 2 · (O_i + S_i) · c / (T_i · (1 + c·T_i)²)
293 ///
294 /// (derived from logit-link binomial Var(S|N) = N p (1−p) and
295 /// d logit(p) / dT = 1/T, scaled by 2 since D = −2 L). Then
296 ///
297 /// I(θ)_{j,k} = Σ_i h_i · J_{i,j} · J_{i,k}.
298 ///
299 /// Returns `None` if the transmission model does not provide an analytical
300 /// Jacobian.
301 pub fn fisher_information(
302 &self,
303 params: &[f64],
304 free_param_indices: &[usize],
305 ) -> Result<Option<FlatMatrix>, FittingError> {
306 let t = self.model.evaluate(params)?;
307 self.validate_inputs(t.len())?;
308 let jac = match self
309 .model
310 .analytical_jacobian(params, free_param_indices, &t)
311 {
312 Some(j) => j,
313 None => return Ok(None),
314 };
315 let n_free = free_param_indices.len();
316 let mut info = FlatMatrix::zeros(n_free, n_free);
317 for (i, ((&t_i, &o_i), &s_i)) in t.iter().zip(self.o.iter()).zip(self.s.iter()).enumerate()
318 {
319 if !self.bin_active(i) {
320 continue;
321 }
322 let h = deviance_curvature(s_i, o_i, t_i, self.c);
323 // Mirror the gradient guard: `deviance_curvature` returns 0
324 // for non-finite `t_i`, but `0.0 * NaN = NaN` would still
325 // poison the Fisher matrix when an FD-built Jacobian has a
326 // NaN cell. Skip the row at h == 0, and skip cells that are
327 // not finite.
328 if h == 0.0 {
329 continue;
330 }
331 for j in 0..n_free {
332 let jij = jac.get(i, j);
333 if !jij.is_finite() {
334 continue;
335 }
336 for k in 0..n_free {
337 let jik = jac.get(i, k);
338 if jik.is_finite() {
339 *info.get_mut(j, k) += h * jij * jik;
340 }
341 }
342 }
343 }
344 Ok(Some(info))
345 }
346
347 /// Finite-difference Fisher information.
348 ///
349 /// Fallback for callers whose transmission model does not implement
350 /// [`FitModel::analytical_jacobian`] — i.e., when
351 /// [`Self::fisher_information`] would return `None`. Builds the
352 /// transmission Jacobian column-by-column via central differences and
353 /// assembles
354 ///
355 /// `I(θ)_{j,k} = Σ_i h_i · J_{i,j} · J_{i,k}`
356 ///
357 /// where `h_i = ∂² D / ∂ T_i²` is the per-bin deviance curvature
358 /// `2·(O_i + S_i)·c / (T_i·(1 + c·T_i)²)` (Fisher-scoring form derived
359 /// from binomial logit-link Var(S | N) = N·p·(1−p) with d logit p / dT
360 /// = 1/T — see the module-level docstring §Model). Returns `Ok(None)`
361 /// only if the base model evaluation itself fails.
362 pub fn fisher_information_fd(
363 &self,
364 params: &mut ParameterSet,
365 fd_step: f64,
366 ) -> Result<Option<FlatMatrix>, FittingError> {
367 let free_idx = params.free_indices();
368 let base_values = params.all_values();
369 let t_base = self.model.evaluate(&base_values)?;
370 self.validate_inputs(t_base.len())?;
371 let n_e = t_base.len();
372 let n_free = free_idx.len();
373 if n_free == 0 {
374 return Ok(Some(FlatMatrix::zeros(0, 0)));
375 }
376 let mut jac = FlatMatrix::zeros(n_e, n_free);
377 for (col, &idx) in free_idx.iter().enumerate() {
378 let original = params.params[idx].value;
379 let step = fd_step * (1.0 + original.abs());
380 params.params[idx].value = original + step;
381 params.params[idx].clamp();
382 let forward_step = params.params[idx].value - original;
383 let t_plus = if forward_step.abs() >= PIVOT_FLOOR {
384 Some(self.model.evaluate(¶ms.all_values())?)
385 } else {
386 None
387 };
388 params.params[idx].value = original - step;
389 params.params[idx].clamp();
390 let backward_step = original - params.params[idx].value;
391 let t_minus = if backward_step.abs() >= PIVOT_FLOOR {
392 Some(self.model.evaluate(¶ms.all_values())?)
393 } else {
394 None
395 };
396 params.params[idx].value = original;
397 let (t_a, t_b, denom) = match (t_plus, t_minus) {
398 (Some(tp), Some(tm)) => (tp, tm, forward_step + backward_step),
399 (Some(tp), None) => (tp, t_base.clone(), forward_step),
400 (None, Some(tm)) => (t_base.clone(), tm, backward_step),
401 (None, None) => continue,
402 };
403 if denom.abs() < PIVOT_FLOOR {
404 continue;
405 }
406 // Per-cell finiteness check. The matching guard in lm.rs
407 // `compute_jacobian` zeroes NaN entries instead of dropping
408 // the column; the same pattern applies here because the
409 // downstream Fisher accumulator below already skips inactive
410 // rows (`bin_active(i)`), so a NaN at a masked / inactive
411 // row must not block the column for active rows. Active-row
412 // NaN is handled by [`deviance_curvature`], which returns 0
413 // on non-finite `t_i` so the assembly stays clean.
414 for i in 0..n_e {
415 let a = t_a[i];
416 let b = t_b[i];
417 if a.is_finite() && b.is_finite() {
418 *jac.get_mut(i, col) = (a - b) / denom;
419 }
420 // else: leave at the zero-default; masked rows are never
421 // read by the active-bin filter in the Fisher loop.
422 }
423 }
424 let mut info = FlatMatrix::zeros(n_free, n_free);
425 for (i, ((&t_i, &o_i), &s_i)) in t_base
426 .iter()
427 .zip(self.o.iter())
428 .zip(self.s.iter())
429 .enumerate()
430 {
431 if !self.bin_active(i) {
432 continue;
433 }
434 let h = deviance_curvature(s_i, o_i, t_i, self.c);
435 // Same guard as the analytical `fisher_information`: avoid
436 // `0.0 * NaN = NaN` poisoning the matrix from NaN Jacobian
437 // cells (per-cell zero default from the FD loop above leaves
438 // most NaN entries as 0, but a stale value from a partial
439 // FD failure must still be defensively skipped).
440 if h == 0.0 {
441 continue;
442 }
443 for j in 0..n_free {
444 let jij = jac.get(i, j);
445 if !jij.is_finite() {
446 continue;
447 }
448 for k in 0..n_free {
449 let jik = jac.get(i, k);
450 if jik.is_finite() {
451 *info.get_mut(j, k) += h * jij * jik;
452 }
453 }
454 }
455 }
456 Ok(Some(info))
457 }
458
459 /// Finite-difference gradient of the deviance.
460 ///
461 /// Central differences on each free parameter. Used as a fallback when
462 /// the model has no analytical Jacobian. `params` is a mutable
463 /// `ParameterSet` so we can respect bounds via `clamp()`.
464 pub fn deviance_gradient_fd(
465 &self,
466 params: &mut ParameterSet,
467 fd_step: f64,
468 ) -> Result<Vec<f64>, FittingError> {
469 let free_idx = params.free_indices();
470 let base_values = params.all_values();
471 let base_d = self.deviance(&base_values)?;
472
473 let mut grad = vec![0.0; free_idx.len()];
474 for (j, &idx) in free_idx.iter().enumerate() {
475 let original = params.params[idx].value;
476 let step = fd_step * (1.0 + original.abs());
477
478 params.params[idx].value = original + step;
479 params.params[idx].clamp();
480 let mut actual_step = params.params[idx].value - original;
481 if actual_step.abs() < PIVOT_FLOOR {
482 // Upper bound blocks forward step: try backward.
483 params.params[idx].value = original - step;
484 params.params[idx].clamp();
485 actual_step = params.params[idx].value - original;
486 if actual_step.abs() < PIVOT_FLOOR {
487 params.params[idx].value = original;
488 continue;
489 }
490 }
491 let perturbed_values = params.all_values();
492 // After the NaN-T contract in `binomial_deviance_term`,
493 // `self.deviance` can legitimately return `Ok(NaN)` when a
494 // probe lands in a region where the model produces a
495 // non-finite transmission. A non-finite `perturbed_d`
496 // divided by `actual_step` would write NaN into `grad[j]`
497 // and poison every subsequent step that consumes the
498 // gradient — symmetric with the `Err` branch below. Treat
499 // both as "this probe is invalid; leave the column at 0".
500 let perturbed_d = match self.deviance(&perturbed_values) {
501 Ok(v) if v.is_finite() => v,
502 _ => {
503 params.params[idx].value = original;
504 continue;
505 }
506 };
507 params.params[idx].value = original;
508 grad[j] = (perturbed_d - base_d) / actual_step;
509 }
510 Ok(grad)
511 }
512}
513
514/// Per-bin binomial deviance term with smooth guards.
515///
516/// Returns `2 · [S·ln(S/(Np)) + O·ln(O/(N(1−p)))]` with the zero-count
517/// convention `x · ln(x / ·) → 0` when `x = 0`.
518///
519/// NaN-T contract (see also [`deviance_weight`] / [`deviance_curvature`]):
520///
521/// - For `0 ≤ T ≤ POISSON_EPSILON` (finite but numerically tiny or zero):
522/// clamps `T` to `POISSON_EPSILON` in the denominator so the optimizer
523/// sees a finite (large) D and a continuous gradient. This is the
524/// "smooth guard" path.
525/// - For **non-finite** `T` (NaN or ±∞): returns `NaN` so the deviance
526/// sum becomes `NaN` and the LM / damped-Fisher trial-step guards
527/// (`Ok(v) if v.is_finite()`) reject the step. This deliberately does
528/// *not* clamp via `f64::max`, because `f64::max(NaN, ε)` returns `ε`
529/// — which would silently masquerade as a valid bin.
530#[inline]
531fn binomial_deviance_term(s: f64, o: f64, t: f64, c: f64) -> f64 {
532 debug_assert!(
533 s.is_finite() && s >= 0.0,
534 "binomial_deviance_term: S must be finite and >= 0, got {s}"
535 );
536 debug_assert!(
537 o.is_finite() && o >= 0.0,
538 "binomial_deviance_term: O must be finite and >= 0, got {o}"
539 );
540 debug_assert!(
541 c.is_finite() && c > 0.0,
542 "binomial_deviance_term: c must be finite and > 0, got {c}"
543 );
544 // `f64::max(NaN, ε)` returns `ε`, so a non-finite T would silently
545 // masquerade as a tiny positive transmission and the deviance would
546 // evaluate to a finite (but meaningless) value that the LM trial-step
547 // guard `Ok(v) if v.is_finite()` would accept. Return NaN so the
548 // deviance sum becomes NaN and the trial step is rejected. The
549 // matching `deviance_weight` / `deviance_curvature` guards return 0,
550 // which keeps the gradient / Fisher accumulators clean rather than
551 // poisoning them with NaN contributions.
552 if !t.is_finite() {
553 return f64::NAN;
554 }
555 let t_safe = t.max(POISSON_EPSILON);
556 let n = s + o;
557 if n <= 0.0 {
558 // Bin has zero counts in both arms — no information, no contribution.
559 return 0.0;
560 }
561 let ct = c * t_safe;
562 // Use a numerically stable form for p. For small cT, p ≈ cT, 1−p ≈ 1.
563 let one_plus_ct = 1.0 + ct;
564 // Expected sample and open-beam counts under profile λ̂.
565 let exp_s = ct / one_plus_ct * n; // = N·p = c·N·T/(1+cT)
566 let exp_o = n / one_plus_ct; // = N·(1−p) = N/(1+cT)
567
568 let term_s = xlogy_ratio(s, exp_s);
569 let term_o = xlogy_ratio(o, exp_o);
570 2.0 * (term_s + term_o)
571}
572
573/// Reject non-finite or negative count arrays at public entry points.
574///
575/// Two distinct failure modes motivate the up-front check:
576///
577/// - **Non-finite (NaN / ±∞).** The per-bin `xlogy_ratio` helper treats
578/// `x <= 0.0` as the zero-count branch and returns 0, but `NaN <= 0.0`
579/// is `false`, so a NaN slips past the branch and propagates
580/// `NaN · ln(NaN / y) = NaN` straight into the deviance sum. The LM
581/// trial-step guard then sees `Ok(NaN)` instead of a clean error.
582/// - **Negative.** `x <= 0.0` *is* true for `x < 0.0`, so the zero-count
583/// branch silently swallows negatives and returns 0 — the deviance
584/// stays finite but the bin is treated as "no data", which is
585/// physically meaningless and conceals the upstream bug (subtraction
586/// artefact in TOF normalisation, signed-int overflow in the loader,
587/// etc.). Negatives never produce NaN, but the "successful" fit
588/// silently discards real data.
589///
590/// Validate up-front so callers get a typed `InvalidConfig` error
591/// pointing at the offending bin instead of either failure mode.
592fn validate_counts(counts: &[f64], field: &'static str) -> Result<(), FittingError> {
593 for (i, &v) in counts.iter().enumerate() {
594 if !v.is_finite() || v < 0.0 {
595 return Err(FittingError::InvalidConfig(format!(
596 "{field}[{i}] must be finite and >= 0, got {v}"
597 )));
598 }
599 }
600 Ok(())
601}
602
603/// `x · ln(x / y)` with the `0 · ln(0 / 0) → 0`, `x · ln(x / 0) → +∞`
604/// conventions. For `y > 0` and `x = 0` the term is 0. For `y = 0` and
605/// `x > 0` we clamp `y` to `POISSON_EPSILON` so the objective stays
606/// finite and continuous.
607#[inline]
608fn xlogy_ratio(x: f64, y: f64) -> f64 {
609 if x <= 0.0 {
610 0.0
611 } else {
612 let y_safe = y.max(POISSON_EPSILON);
613 x * (x / y_safe).ln()
614 }
615}
616
617/// Per-bin ∂D/∂T.
618///
619/// ∂D/∂T = −2 · (S − O·c·T) / (T · (1 + c·T))
620///
621/// When `T ≤ ε`, uses a linear extrapolation from `T = ε` so the gradient
622/// stays finite and continuous across the boundary (matching the clamping
623/// done in [`binomial_deviance_term`]).
624#[inline]
625fn deviance_weight(s: f64, o: f64, t: f64, c: f64) -> f64 {
626 // A non-finite T must not be folded into the gradient accumulator.
627 // `f64::max(NaN, ε)` returns `ε`, which would turn a NaN bin into a
628 // finite gradient contribution scaled by the Jacobian and silently
629 // steer the optimizer. Skip the bin (return 0) — the matching
630 // `binomial_deviance_term` returns NaN so the step is rejected by
631 // the trial-guard, but the gradient stays clean in case the caller
632 // is using it for diagnostics on a partially-bad grid.
633 if !t.is_finite() {
634 return 0.0;
635 }
636 let t_safe = t.max(POISSON_EPSILON);
637 let one_plus_ct = 1.0 + c * t_safe;
638 -2.0 * (s - o * c * t_safe) / (t_safe * one_plus_ct)
639}
640
641/// Per-bin ∂²D/∂T² using the expected-info (Fisher) form.
642///
643/// Under the model, Var(S | N) = N · p · (1 − p) = N · cT / (1+cT)². With
644/// d logit(p) / dT = 1/T, the Fisher info on T is
645///
646/// I_TT = N · c / (T · (1 + c·T)²)
647///
648/// and ∂²D/∂T² = 2 · I_TT (since D = −2 · L_c).
649#[inline]
650fn deviance_curvature(s: f64, o: f64, t: f64, c: f64) -> f64 {
651 // See the matching guard in [`deviance_weight`]. A non-finite T
652 // would otherwise contribute a huge spurious curvature via
653 // `f64::max(NaN, ε) -> ε`, inflating the diagonal of the Fisher
654 // matrix and underestimating the corresponding parameter
655 // uncertainty (covariance = I⁻¹ entries shrink as I grows).
656 if !t.is_finite() {
657 return 0.0;
658 }
659 let t_safe = t.max(POISSON_EPSILON);
660 let n = s + o;
661 let one_plus_ct = 1.0 + c * t_safe;
662 2.0 * n * c / (t_safe * one_plus_ct * one_plus_ct)
663}
664
665// ======================================================================
666// joint_poisson_fit — two-stage solver (damped Fisher + Nelder-Mead polish)
667// ======================================================================
668
669use crate::lm::{invert_matrix, solve_damped_system};
670use crate::nelder_mead::{NelderMeadConfig, nelder_mead_minimize};
671
672/// Configuration for [`joint_poisson_fit`].
673#[derive(Debug, Clone)]
674pub struct JointPoissonFitConfig {
675 /// Maximum number of damped-Fisher iterations in stage 1.
676 pub max_iter: usize,
677 /// Initial damping factor (Marquardt λ) on the Fisher matrix diagonal.
678 pub lambda_init: f64,
679 /// Multiplicative factor to increase λ on a rejected step.
680 pub lambda_up: f64,
681 /// Multiplicative factor to decrease λ on an accepted step.
682 pub lambda_down: f64,
683 /// Armijo sufficient-decrease coefficient.
684 pub armijo_c: f64,
685 /// Backtracking factor during line search.
686 pub backtrack: f64,
687 /// Convergence tolerance on relative deviance change.
688 pub tol_d: f64,
689 /// Convergence tolerance on normalized parameter step.
690 pub tol_param: f64,
691 /// Finite-difference step for gradient fallback.
692 pub fd_step: f64,
693 /// Enable Nelder-Mead polish after stage 1.
694 ///
695 /// Default `false` as of #486. The polish tolerances
696 /// (`xatol = 1e-9, fatol = 1e-10`) were originally matched to the
697 /// EG5 synthetic benchmark (memo 35 §P2.1) where D stays O(1), so
698 /// `fatol` is physically meaningful. On real-data regimes where
699 /// D saturates at 10⁴–10⁵ (un-modelled upstream physics —
700 /// memo 35 §P3/§P4), `fatol / D` drops below f64 ULP and polish
701 /// cannot self-terminate — it burns its full `max_iter = 5000`
702 /// every fit at 70–260× wall cost, and the three-scenario
703 /// ablation on real VENUS Hf 120-min data (issue #486) showed
704 /// the resulting parameter shift is ≤ 0.35 Fisher σ on every
705 /// parameter in every scenario — i.e. below the solver's own
706 /// reported uncertainty floor.
707 ///
708 /// The polish mechanism itself is sound (self-terminates cleanly
709 /// on synthetic D≈1 data per ablation S3); only the absolute
710 /// tolerance defaults are mis-calibrated for real counts data.
711 /// A future scale-aware rescale (`fatol_rel` vs `D_stage1`) can
712 /// re-enable polish as a useful opt-in refinement.
713 ///
714 /// Set this to `true` (via `with_counts_enable_polish(Some(true))`
715 /// at the pipeline level) when you specifically want the polish
716 /// stage on a synthetic / clean-data scenario where the absolute
717 /// tolerance defaults are physically meaningful.
718 pub enable_polish: bool,
719 /// Polish (Nelder-Mead) configuration. Used only when
720 /// `enable_polish == true`. Default `xatol = 1e-9`, `fatol = 1e-10`
721 /// match the EG5 synthetic benchmark tolerances from memo 35 §P2.1
722 /// — physically meaningful when `D ≈ 1` (clean data) but sub-f64-
723 /// ULP on real counts where `D ≈ 10⁴`–`10⁵`, which is why
724 /// `enable_polish` defaults to `false`. See #486.
725 pub polish: NelderMeadConfig,
726 /// Compute and return the Fisher covariance and parameter uncertainties.
727 pub compute_covariance: bool,
728}
729
730impl Default for JointPoissonFitConfig {
731 fn default() -> Self {
732 Self {
733 max_iter: 200,
734 lambda_init: 1e-3,
735 lambda_up: 10.0,
736 lambda_down: 0.1,
737 armijo_c: 1e-4,
738 backtrack: 0.5,
739 tol_d: 1e-8,
740 tol_param: 1e-8,
741 fd_step: 1e-6,
742 // #486: flipped from `true` to `false` after a three-scenario
743 // ablation on real VENUS data showed polish burning full
744 // `max_iter = 5000` at 70-260× wall cost for ≤ 0.35 Fisher σ
745 // parameter movement. The absolute tolerances below are
746 // physically meaningful for synthetic (D ≈ 1) benchmarks and
747 // dead on real counts data (D ≈ 10⁵). Opt in via
748 // `UnifiedFitConfig::with_counts_enable_polish(Some(true))`
749 // when you specifically want the polish stage. See the
750 // field doc on `enable_polish` for details.
751 enable_polish: false,
752 polish: NelderMeadConfig {
753 // Tolerances tuned for the EG5 synthetic regime (memo 35
754 // §P2.1) — `fatol = 1e-10` vs D ≈ 1 is a physically
755 // meaningful "deviance isn't budging" check. On real
756 // counts data where D ≈ 10⁵ the same absolute value is
757 // sub-ULP; polish can't self-terminate and is disabled
758 // by the default above. A future scale-aware rescale
759 // (`fatol_rel` vs D_stage1) is tracked as a follow-up.
760 xatol: 1e-9,
761 fatol: 1e-10,
762 max_iter: 5000,
763 initial_step_frac: 0.02,
764 initial_step_abs: 1e-4,
765 },
766 compute_covariance: true,
767 }
768 }
769}
770
771/// Outcome of [`joint_poisson_fit`].
772#[derive(Debug, Clone)]
773pub struct JointPoissonResult {
774 /// Final deviance D at the fitted parameters.
775 pub deviance: f64,
776 /// D / (n − k). Primary GOF statistic per memo 35 §P1.2.
777 pub deviance_per_dof: f64,
778 /// Number of data bins on the configured grid (n). This is the
779 /// total bin count; when a fit-energy-range mask is in effect, the
780 /// count of bins that actually contributed to the cost function is
781 /// reported separately as [`Self::n_active`].
782 pub n_data: usize,
783 /// Number of *active* data bins — equal to `n_data` when no mask is
784 /// set, or the count of `true` entries in the objective's
785 /// `active_mask` otherwise. The deviance / dof ratio uses
786 /// `(n_active − n_free)` so reduced deviance is unbiased when a
787 /// fit-energy-range mask is in effect (SAMMY EMIN/EMAX semantics, #514).
788 pub n_active: usize,
789 /// Number of free parameters (k).
790 pub n_free: usize,
791 /// Iterations performed in the damped-Fisher stage.
792 pub gn_iterations: usize,
793 /// Iterations performed by the Nelder-Mead polish stage (0 if disabled).
794 pub polish_iterations: usize,
795 /// `true` when the stage-1 (damped Fisher) optimizer met its `tol_d`
796 /// and `tol_param` criteria before hitting `max_iter`.
797 pub gn_converged: bool,
798 /// `true` when the Nelder-Mead polish met `xatol` and `fatol` before
799 /// `max_iter` (always `false` if `enable_polish == false`).
800 pub polish_converged: bool,
801 /// `true` when the polish step lowered the deviance below the stage-1
802 /// best value. Useful diagnostic — if polish improved D materially,
803 /// stage 1 likely stalled.
804 pub polish_improved: bool,
805 /// Final parameter values (all parameters, including fixed).
806 pub params: Vec<f64>,
807 /// Inverse Fisher covariance of free parameters (n_free × n_free),
808 /// computed at the final θ. `None` if the Fisher matrix was singular
809 /// or `compute_covariance == false`.
810 pub covariance: Option<FlatMatrix>,
811 /// `√diag(covariance)` for each free parameter, in free-index order.
812 pub uncertainties: Option<Vec<f64>>,
813}
814
815/// Two-stage joint-Poisson fit: damped Fisher stage followed by
816/// Nelder-Mead polish.
817///
818/// **Memo 35 §P1 + §P2 requirements** this function satisfies:
819///
820/// - Minimizes the **conditional binomial deviance** `D(θ)`
821/// ([`JointPoissonObjective::deviance`]), not fixed-flux Poisson NLL.
822/// - Reports `D / (n − k)` as the primary GOF (P1.2).
823/// - Honours an **explicit `c = Q_s/Q_ob`** stored in the objective (P1.3).
824/// - Runs Nelder-Mead **polish** after the gradient stage to escape the
825/// EG2-S1 C_full initial-point stall (P2.1).
826/// - Exposes `gn_converged` and `polish_converged` separately so callers
827/// do not rely on a single "success" flag — acceptance is meant to come
828/// from the deviance value (P2.3).
829///
830/// The damped-Fisher stage uses LM-style acceptance: a step is accepted if
831/// it satisfies an Armijo condition on D; on rejection, λ is increased and
832/// the step is recomputed. Bounds are enforced via projection (clamp).
833pub fn joint_poisson_fit(
834 objective: &JointPoissonObjective<'_>,
835 params: &mut ParameterSet,
836 config: &JointPoissonFitConfig,
837) -> Result<JointPoissonResult, FittingError> {
838 let n_data = objective.n_data();
839 if n_data == 0 {
840 return Err(FittingError::EmptyData);
841 }
842
843 // Validate `o` / `s` length and `c` up-front at the public entry
844 // point. The inner per-bin helpers (`binomial_deviance_term`,
845 // `deviance_from_transmission`) use `debug_assert!` only, which is a
846 // no-op in release builds. Without these hard checks:
847 // - A length mismatch in `o` vs `s` silently truncates via `.zip()`,
848 // minimising deviance on a sub-range of bins.
849 // - A non-positive or non-finite `c` produces finite garbage
850 // (e.g. zero `cT`, NaN denominators) that the LM happily descends.
851 // - A NaN / negative `o[i]` or `s[i]` would slip past the inner
852 // `xlogy_ratio` zero-clamp (`x <= 0.0` swallows negatives, but
853 // `NaN <= 0.0` is `false` so a NaN count bleeds straight into the
854 // log and out into the deviance sum).
855 // All surface as "the fit converged" with bogus parameter values —
856 // exactly the failure mode the trial-step guard cannot catch because
857 // the deviance value is finite.
858 if objective.s.len() != n_data {
859 return Err(FittingError::LengthMismatch {
860 expected: n_data,
861 actual: objective.s.len(),
862 field: "sample_counts",
863 });
864 }
865 if !objective.c.is_finite() || objective.c <= 0.0 {
866 return Err(FittingError::InvalidConfig(format!(
867 "proton-charge ratio c = Q_s/Q_ob must be finite and > 0, got {}",
868 objective.c
869 )));
870 }
871 validate_counts(objective.o, "open_beam_counts")?;
872 validate_counts(objective.s, "sample_counts")?;
873
874 // Validate active-mask length up-front, mirroring the LM solver's
875 // length-mismatch early-return (#514). A debug-assert deep in the
876 // deviance routines would silently pass through in release builds
877 // with a length mismatch, causing out-of-bounds index reads when
878 // the masked accumulator iterates `o`/`s`/`mask` together.
879 if let Some(m) = objective.active_mask
880 && m.len() != n_data
881 {
882 return Err(FittingError::LengthMismatch {
883 expected: n_data,
884 actual: m.len(),
885 field: "active_mask",
886 });
887 }
888
889 // SAMMY EMIN/EMAX-equivalent fit-energy-range (#514): zero active
890 // bins means the user's `[E_min, E_max]` does not overlap the
891 // configured grid. No data contributes to the deviance — return
892 // non-converged with NaN before falling through. Without this
893 // the all-bins-masked path would compute deviance = 0 (sum over
894 // zero rows) which combined with `n_free == 0` (all-fixed params)
895 // would report `gn_converged: true, deviance: 0` from a fit that
896 // saw no data.
897 let n_free_initial = params.n_free();
898 let n_active_initial = objective.n_active();
899 if n_active_initial == 0 {
900 return Ok(JointPoissonResult {
901 deviance: f64::NAN,
902 deviance_per_dof: f64::NAN,
903 n_data,
904 n_active: 0,
905 n_free: n_free_initial,
906 gn_iterations: 0,
907 polish_iterations: 0,
908 gn_converged: false,
909 polish_converged: false,
910 polish_improved: false,
911 params: params.all_values(),
912 covariance: None,
913 uncertainties: None,
914 });
915 }
916
917 // Underdetermined-check: when a fit-energy-range mask leaves fewer
918 // active bins than free parameters, the problem is rank-deficient
919 // and any deviance / dof ratio would be deceptive (the previous
920 // `.max(1)` divisor produced a finite-looking deviance-per-dof for
921 // empty / too-narrow masks). Mirror LM's behaviour at
922 // `lm.rs:578-588`: return a non-converged result up-front, before
923 // wasting cycles on the damped-Fisher stage.
924 if n_active_initial < n_free_initial {
925 return Ok(JointPoissonResult {
926 deviance: f64::NAN,
927 deviance_per_dof: f64::NAN,
928 n_data,
929 n_active: n_active_initial,
930 n_free: n_free_initial,
931 gn_iterations: 0,
932 polish_iterations: 0,
933 gn_converged: false,
934 polish_converged: false,
935 polish_improved: false,
936 params: params.all_values(),
937 covariance: None,
938 uncertainties: None,
939 });
940 }
941
942 // Stage 1: damped Fisher with Armijo backtracking.
943 let stage1 = damped_fisher_stage(objective, params, config)?;
944
945 // Capture stage-1 best.
946 let best_d_stage1 = stage1.deviance;
947 let gn_iterations = stage1.iterations;
948 let gn_converged = stage1.converged;
949
950 // Stage 2: Nelder-Mead polish on free parameters, seeded from stage-1 θ.
951 //
952 // Guard against the all-fixed configuration: `nelder_mead_minimize`
953 // requires a non-empty `x0` (asserts in `nelder_mead.rs`). When every
954 // parameter is fixed there is nothing to polish, so skip stage 2 and
955 // leave the polish flags at their default `false` values. This path
956 // is reachable from pipeline callers that pin all params and set
957 // `with_counts_enable_polish(Some(true))`.
958 //
959 // Also short-circuit polish when stage 1 ended on a non-finite
960 // deviance: there is no meaningful starting deviance to refine, and
961 // the acceptance test `nm.fun < best_d_stage1` would degrade to
962 // `finite < NaN == false` (discarding the NM result) while
963 // `nm.self_converged` could still be `true`, leaking a spurious
964 // converged flag together with a NaN final deviance. Mirrors the
965 // LM `n_free == 0` early-return at `lm.rs:584-607`, which refuses to
966 // report a converged fit when the model emits NaN at active bins.
967 let mut polish_iterations = 0usize;
968 let mut polish_converged = false;
969 let mut polish_improved = false;
970 let free_idx = params.free_indices();
971 if config.enable_polish && !free_idx.is_empty() && best_d_stage1.is_finite() {
972 let bounds: Vec<(f64, f64)> = free_idx
973 .iter()
974 .map(|&i| (params.params[i].lower, params.params[i].upper))
975 .collect();
976 let x0: Vec<f64> = free_idx.iter().map(|&i| params.params[i].value).collect();
977
978 // Snapshot fixed parameters so the closure can rebuild the full
979 // parameter vector for each evaluation.
980 let all_values_snapshot = params.all_values();
981
982 let obj_closure = |x: &[f64]| -> Result<f64, FittingError> {
983 let mut all = all_values_snapshot.clone();
984 for (j, &idx) in free_idx.iter().enumerate() {
985 all[idx] = x[j];
986 }
987 objective.deviance(&all)
988 };
989 let nm = nelder_mead_minimize(obj_closure, &x0, Some(&bounds), &config.polish)?;
990 polish_iterations = nm.iterations;
991 polish_converged = nm.self_converged;
992 if nm.fun < best_d_stage1 {
993 polish_improved = true;
994 // Commit polish result to the parameter set.
995 for (j, &idx) in free_idx.iter().enumerate() {
996 params.params[idx].value = nm.x[j];
997 params.params[idx].clamp();
998 }
999 }
1000 }
1001
1002 let final_values = params.all_values();
1003 let final_deviance = objective.deviance(&final_values)?;
1004 let n_free = params.n_free();
1005 // Active-bin masking (SAMMY EMIN/EMAX): when a fit-energy-range mask
1006 // is in effect, dof must use the count of bins that contributed to
1007 // the deviance — otherwise deviance-per-dof is biased low by the
1008 // ratio (n_active / n_data). The `n_active < n_free` case has
1009 // already been short-circuited above; here `n_active >= n_free`,
1010 // so `dof` is non-negative and exactly-determined fits
1011 // (`n_active == n_free`) report `deviance_per_dof = NaN` (0/0)
1012 // as in LM (`lm.rs:784`).
1013 let n_active = objective.n_active();
1014 let dof = n_active.saturating_sub(n_free);
1015 let deviance_per_dof = if dof > 0 {
1016 final_deviance / dof as f64
1017 } else {
1018 f64::NAN
1019 };
1020
1021 // Covariance from inverse Fisher at the final θ. Uses the analytical
1022 // Jacobian when the transmission model provides one; otherwise falls
1023 // back to finite-difference Jacobian assembled into the deviance-
1024 // Hessian form — so callers always get uncertainties for identifiable
1025 // parameters.
1026 //
1027 // **Scale note (covariance vs Newton step).** `fisher_information`
1028 // assembles `H_D = Σ h_i · J·J^T` with `h_i = ∂² D / ∂ T_i² = 2 · I_TT_i`
1029 // (see [`deviance_curvature`]). This `2·I` form is exactly what the
1030 // damped-Fisher Newton step needs, since stepping on D with
1031 // `Δθ = -H_D^{-1} · ∇D = -(2I)^{-1} · (-2 ∇L) = I^{-1} · ∇L`
1032 // recovers the Fisher-scoring direction on the log-likelihood L.
1033 //
1034 // For the asymptotic MLE covariance, however, the Cramér-Rao bound is
1035 // `Cov(θ̂) = I^{-1}`, NOT `H_D^{-1} = (2I)^{-1} = I^{-1}/2`. Inverting
1036 // `H_D` and using it directly would under-report variance by 2× and
1037 // standard errors by √2 × — a real bug caught in review. We rescale
1038 // the inverse here: `I^{-1} = 2 · H_D^{-1}`.
1039 let (covariance, uncertainties) = if config.compute_covariance {
1040 let free_idx = params.free_indices();
1041 let info_opt = match objective.fisher_information(&final_values, &free_idx)? {
1042 Some(info) => Some(info),
1043 None => objective.fisher_information_fd(params, config.fd_step)?,
1044 };
1045 match info_opt {
1046 Some(info) => match invert_matrix(&info) {
1047 Some(mut cov) => {
1048 // Rescale: invert_matrix returned (2I)^{-1}; multiply
1049 // every entry by 2 to obtain I^{-1}.
1050 for v in cov.data.iter_mut() {
1051 *v *= 2.0;
1052 }
1053 let u: Vec<f64> = (0..cov.nrows)
1054 .map(|i| {
1055 let v = cov.get(i, i);
1056 if v > 0.0 { v.sqrt() } else { f64::NAN }
1057 })
1058 .collect();
1059 (Some(cov), Some(u))
1060 }
1061 None => (None, None),
1062 },
1063 None => (None, None),
1064 }
1065 } else {
1066 (None, None)
1067 };
1068
1069 Ok(JointPoissonResult {
1070 deviance: final_deviance,
1071 deviance_per_dof,
1072 n_data,
1073 n_active,
1074 n_free,
1075 gn_iterations,
1076 polish_iterations,
1077 gn_converged,
1078 polish_converged,
1079 polish_improved,
1080 params: final_values,
1081 covariance,
1082 uncertainties,
1083 })
1084}
1085
1086/// Stage 1 output.
1087struct Stage1Output {
1088 deviance: f64,
1089 iterations: usize,
1090 converged: bool,
1091}
1092
1093/// Damped-Fisher stage (Gauss-Newton / Marquardt on the deviance).
1094///
1095/// Mirrors the structure of `lm.rs` but on the joint-Poisson objective.
1096/// Falls back to finite-difference gradient when the model has no
1097/// analytical Jacobian.
1098fn damped_fisher_stage(
1099 objective: &JointPoissonObjective<'_>,
1100 params: &mut ParameterSet,
1101 config: &JointPoissonFitConfig,
1102) -> Result<Stage1Output, FittingError> {
1103 let mut lambda = config.lambda_init;
1104 let mut iter = 0usize;
1105 let mut converged = false;
1106
1107 let mut all_vals = params.all_values();
1108 let mut d_current = objective.deviance(&all_vals)?;
1109
1110 while iter < config.max_iter {
1111 iter += 1;
1112 let free_idx = params.free_indices();
1113 let n_free = free_idx.len();
1114 if n_free == 0 {
1115 // All parameters fixed: we are not optimizing; convergence is
1116 // well-defined only if the already-computed deviance at the
1117 // current parameters is finite. If the model returned
1118 // non-finite transmission, `binomial_deviance_term` propagates
1119 // that as NaN deviance (see the non-finite-T contract documented
1120 // on `binomial_deviance_term`), and a non-finite deviance cannot
1121 // be reported as a converged fit. LM applies the same guard in
1122 // the `n_free == 0` branch of `levenberg_marquardt_with_mask`;
1123 // the matching LM regression is `test_all_fixed_params_nan_model`.
1124 converged = d_current.is_finite();
1125 break;
1126 }
1127
1128 // Gradient (analytical if available, FD otherwise).
1129 let grad = match objective.deviance_gradient_analytical(&all_vals, &free_idx)? {
1130 Some(g) => g,
1131 None => objective.deviance_gradient_fd(params, config.fd_step)?,
1132 };
1133 // Fisher information (Gauss-Newton curvature). If absent, use a
1134 // diagonal identity fallback scaled by gradient magnitude — this
1135 // degenerates the stage into projected gradient descent, which is
1136 // exactly how `poisson.rs` behaves in the FD regime.
1137 let info = match objective.fisher_information(&all_vals, &free_idx)? {
1138 Some(m) => m,
1139 None => {
1140 let mut ident = FlatMatrix::zeros(n_free, n_free);
1141 for i in 0..n_free {
1142 *ident.get_mut(i, i) = 1.0;
1143 }
1144 ident
1145 }
1146 };
1147 // Solve (I + λ diag(I)) δ = -g.
1148 let neg_grad: Vec<f64> = grad.iter().map(|&g| -g).collect();
1149 let step = match solve_damped_system(&info, &neg_grad, lambda) {
1150 Some(s) => s,
1151 None => {
1152 // Singular Fisher at current θ. Increase damping and retry
1153 // on the next iteration.
1154 lambda *= config.lambda_up;
1155 if lambda > 1e16 {
1156 break;
1157 }
1158 continue;
1159 }
1160 };
1161
1162 // Armijo line search with projection.
1163 let grad_dot_step = grad
1164 .iter()
1165 .zip(step.iter())
1166 .map(|(&g, &s)| g * s)
1167 .sum::<f64>();
1168 // If the step isn't a descent direction w.r.t. D, flip sign (fallback
1169 // to negative gradient direction).
1170 let effective_step: Vec<f64> = if grad_dot_step >= 0.0 {
1171 grad.iter().map(|&g| -g).collect()
1172 } else {
1173 step
1174 };
1175
1176 let mut alpha = 1.0;
1177 let mut accepted = false;
1178 let d0 = d_current;
1179 let mut trial_vals = all_vals.clone();
1180 for _ in 0..50 {
1181 for (j, &idx) in free_idx.iter().enumerate() {
1182 trial_vals[idx] = all_vals[idx] + alpha * effective_step[j];
1183 }
1184 // Project onto bounds.
1185 for &idx in free_idx.iter() {
1186 let lo = params.params[idx].lower;
1187 let hi = params.params[idx].upper;
1188 if trial_vals[idx] < lo {
1189 trial_vals[idx] = lo;
1190 }
1191 if trial_vals[idx] > hi {
1192 trial_vals[idx] = hi;
1193 }
1194 }
1195 let d_trial = match objective.deviance(&trial_vals) {
1196 Ok(v) if v.is_finite() => v,
1197 _ => f64::INFINITY,
1198 };
1199 // Armijo condition: f(x+αp) ≤ f(x) + c·α·⟨g, p⟩ (descent). When
1200 // we flipped to -grad above, ⟨g, p⟩ = -||g||² < 0.
1201 let gdotp = grad
1202 .iter()
1203 .zip(effective_step.iter())
1204 .map(|(&g, &s)| g * s)
1205 .sum::<f64>();
1206 if d_trial <= d0 + config.armijo_c * alpha * gdotp {
1207 accepted = true;
1208 break;
1209 }
1210 alpha *= config.backtrack;
1211 if alpha < 1e-16 {
1212 break;
1213 }
1214 }
1215
1216 if accepted {
1217 // Commit step.
1218 for &idx in free_idx.iter() {
1219 params.params[idx].value = trial_vals[idx];
1220 params.params[idx].clamp();
1221 }
1222 let rel_change =
1223 (d_current - objective.deviance(&trial_vals)?) / d_current.abs().max(1.0);
1224 all_vals = params.all_values();
1225 let new_d = objective.deviance(&all_vals)?;
1226 let step_norm_sq = effective_step
1227 .iter()
1228 .map(|&s| (alpha * s).powi(2))
1229 .sum::<f64>();
1230 let step_norm = step_norm_sq.sqrt();
1231 d_current = new_d;
1232 lambda = (lambda * config.lambda_down).max(1e-16);
1233
1234 if rel_change.abs() < config.tol_d && step_norm < config.tol_param {
1235 converged = true;
1236 break;
1237 }
1238 } else {
1239 // Rejected: increase damping and try again.
1240 lambda *= config.lambda_up;
1241 if lambda > 1e16 {
1242 break;
1243 }
1244 }
1245 }
1246
1247 Ok(Stage1Output {
1248 deviance: d_current,
1249 iterations: iter,
1250 converged,
1251 })
1252}
1253
1254#[cfg(test)]
1255mod tests {
1256 use super::*;
1257 use crate::parameters::FitParameter;
1258
1259 // ------------------------------------------------------------------
1260 // Test fixtures
1261 // ------------------------------------------------------------------
1262
1263 /// A constant-transmission model: T_i = θ_0 for all i. Useful for
1264 /// testing the profile λ̂ formula and deviance / gradient in isolation.
1265 struct ConstModel {
1266 n_e: usize,
1267 }
1268
1269 impl FitModel for ConstModel {
1270 fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1271 Ok(vec![params[0]; self.n_e])
1272 }
1273
1274 fn analytical_jacobian(
1275 &self,
1276 _params: &[f64],
1277 free_param_indices: &[usize],
1278 y_current: &[f64],
1279 ) -> Option<FlatMatrix> {
1280 let n_e = y_current.len();
1281 let n_free = free_param_indices.len();
1282 let mut jac = FlatMatrix::zeros(n_e, n_free);
1283 // ∂T/∂θ_0 = 1 for all i, and 0 for any other parameter.
1284 for i in 0..n_e {
1285 for (j, &pi) in free_param_indices.iter().enumerate() {
1286 *jac.get_mut(i, j) = if pi == 0 { 1.0 } else { 0.0 };
1287 }
1288 }
1289 Some(jac)
1290 }
1291 }
1292
1293 /// A linear-in-E model: T_i = θ_0 − θ_1 · e_i (Beer-Lambert surrogate).
1294 /// Used for the analytical-vs-FD gradient check and profile tests with
1295 /// non-trivial Jacobian.
1296 struct LinearModel<'a> {
1297 e: &'a [f64],
1298 }
1299
1300 impl<'a> FitModel for LinearModel<'a> {
1301 fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1302 Ok(self
1303 .e
1304 .iter()
1305 .map(|&ei| (params[0] - params[1] * ei).max(POISSON_EPSILON))
1306 .collect())
1307 }
1308
1309 fn analytical_jacobian(
1310 &self,
1311 _params: &[f64],
1312 free_param_indices: &[usize],
1313 y_current: &[f64],
1314 ) -> Option<FlatMatrix> {
1315 let n_e = y_current.len();
1316 let n_free = free_param_indices.len();
1317 let mut jac = FlatMatrix::zeros(n_e, n_free);
1318 for i in 0..n_e {
1319 for (j, &pi) in free_param_indices.iter().enumerate() {
1320 *jac.get_mut(i, j) = match pi {
1321 0 => 1.0,
1322 1 => -self.e[i],
1323 _ => 0.0,
1324 };
1325 }
1326 }
1327 Some(jac)
1328 }
1329 }
1330
1331 // ------------------------------------------------------------------
1332 // (a) Profile λ̂ closed form matches the score-equation bisection root.
1333 // ------------------------------------------------------------------
1334 #[test]
1335 fn test_profile_lambda_closed_form_matches_bisection() {
1336 // For each bin independently, score(λ) = (O+S)/λ − (1/c + T) = 0
1337 // has the unique positive root λ̂ = c(O+S)/(1+cT). Bisect on
1338 // [1e-10, 1e12] and verify agreement to 1e-9.
1339 let cases = [
1340 (50.0_f64, 5.0_f64, 0.5_f64, 1.0_f64),
1341 (1000.0, 900.0, 0.9, 5.98),
1342 (10.0, 1.0, 0.1, 2.0),
1343 (0.0, 5.0, 0.25, 1.5), // O=0 edge
1344 (5.0, 0.0, 0.75, 3.0), // S=0 edge
1345 ];
1346 for (o, s, t, c) in cases {
1347 let model = ConstModel { n_e: 1 };
1348 let obj = JointPoissonObjective {
1349 model: &model,
1350 o: &[o],
1351 s: &[s],
1352 c,
1353 active_mask: None,
1354 };
1355 let closed = obj.profile_lambda(t, o, s);
1356
1357 // Bisection root of score(λ) = (O+S)/λ − (1/c + T).
1358 let score = |lam: f64| (o + s) / lam - (1.0 / c + t);
1359 let (mut lo, mut hi) = (1e-10, 1e12);
1360 // score is monotonically decreasing in λ, score(lo) > 0, score(hi) < 0.
1361 assert!(score(lo) >= 0.0);
1362 assert!(score(hi) <= 0.0);
1363 for _ in 0..200 {
1364 let mid = 0.5 * (lo + hi);
1365 if score(mid) > 0.0 {
1366 lo = mid;
1367 } else {
1368 hi = mid;
1369 }
1370 }
1371 let bisect = 0.5 * (lo + hi);
1372 let rel_err = ((closed - bisect) / bisect).abs();
1373 assert!(
1374 rel_err < 1e-9,
1375 "profile λ̂ mismatch: closed={closed} bisect={bisect} rel_err={rel_err}"
1376 );
1377 }
1378 }
1379
1380 // ------------------------------------------------------------------
1381 // (b) D = 0 at exact match of expected counts.
1382 // ------------------------------------------------------------------
1383 #[test]
1384 fn test_deviance_zero_at_exact_match() {
1385 // Construct a model where S_i = λ·T_i, O_i = λ/c exactly for integer
1386 // choices, then verify D < 1e-8. With T=0.5, c=2, λ=200: S=100,
1387 // O=100 per bin; p = 2*0.5/(1+1) = 0.5; Np = (O+S)/2 = 100 = S;
1388 // N(1-p) = 100 = O, so both logs are zero and D = 0.
1389 let t_val = 0.5;
1390 let c = 2.0;
1391 let n_bins = 5;
1392 let o = vec![100.0; n_bins];
1393 let s = vec![100.0; n_bins];
1394 let t = vec![t_val; n_bins];
1395 let model = ConstModel { n_e: n_bins };
1396 let obj = JointPoissonObjective {
1397 model: &model,
1398 o: &o,
1399 s: &s,
1400 c,
1401 active_mask: None,
1402 };
1403 let d = obj.deviance_from_transmission(&t).unwrap();
1404 assert!(d.abs() < 1e-8, "D should be ≈ 0 at exact match, got {d}");
1405
1406 // Also verify via parameter evaluation (model returns constant T).
1407 let d_via_params = obj.deviance(&[t_val]).unwrap();
1408 assert!(d_via_params.abs() < 1e-8);
1409 }
1410
1411 // ------------------------------------------------------------------
1412 // (c) Analytical gradient matches finite-difference.
1413 // ------------------------------------------------------------------
1414 #[test]
1415 fn test_deviance_gradient_matches_fd() {
1416 // Use the linear model T = θ_0 − θ_1 · E with noise-free synthetic
1417 // counts. Compute analytical gradient via chain rule and FD
1418 // gradient via re-evaluation; they must agree.
1419 let e: Vec<f64> = (0..20).map(|i| 0.1 + 0.05 * i as f64).collect();
1420 let theta_true = [0.95_f64, 0.1_f64];
1421 let c = 3.0;
1422 let lam = 500.0;
1423
1424 // Generate noise-free expected counts.
1425 let model = LinearModel { e: &e };
1426 let t_true = model.evaluate(&theta_true).unwrap();
1427 let o: Vec<f64> = t_true.iter().map(|_| lam / c).collect();
1428 let s: Vec<f64> = t_true.iter().map(|&ti| lam * ti).collect();
1429
1430 let obj = JointPoissonObjective {
1431 model: &model,
1432 o: &o,
1433 s: &s,
1434 c,
1435 active_mask: None,
1436 };
1437
1438 // Evaluate gradient at a point slightly off truth so it is nonzero.
1439 let theta_eval = [0.80_f64, 0.15_f64];
1440 let free_idx = vec![0, 1];
1441
1442 let g_analytical = obj
1443 .deviance_gradient_analytical(&theta_eval, &free_idx)
1444 .unwrap()
1445 .expect("LinearModel provides analytical jacobian");
1446
1447 // Central-difference gradient.
1448 let eps = 1e-6;
1449 let mut g_fd = [0.0_f64; 2];
1450 for j in 0..2 {
1451 let mut tp = theta_eval;
1452 let mut tm = theta_eval;
1453 tp[j] += eps;
1454 tm[j] -= eps;
1455 let dp = obj.deviance(&tp).unwrap();
1456 let dm = obj.deviance(&tm).unwrap();
1457 g_fd[j] = (dp - dm) / (2.0 * eps);
1458 }
1459
1460 for (a, f) in g_analytical.iter().zip(g_fd.iter()) {
1461 let rel = ((a - f) / f.abs().max(1e-6)).abs();
1462 assert!(
1463 rel < 1e-4,
1464 "analytical vs FD gradient disagree: analytical={a} fd={f} rel={rel}"
1465 );
1466 }
1467 }
1468
1469 // ------------------------------------------------------------------
1470 // (d) D/(n-k) asymptote on synthetic joint-Poisson data at matched
1471 // model — single free parameter θ_0 = T, use 1D grid search to
1472 // recover it, verify D/(n-1) ≈ 1 and density bias < 1%.
1473 // ------------------------------------------------------------------
1474 #[test]
1475 fn test_deviance_per_dof_asymptote() {
1476 // Deterministic generator (xorshift) so the test is reproducible.
1477 // `Xorshift` is defined at the module level below — Rust item order
1478 // is not significant inside a module.
1479 let n_bins = 200;
1480 let t_true = 0.35_f64;
1481 let c = 2.0;
1482 let lam = 50.0;
1483 let n_reps = 30;
1484
1485 let mut d_per_dof_samples = Vec::with_capacity(n_reps);
1486 let mut bias_samples = Vec::with_capacity(n_reps);
1487 let mut rng = Xorshift(0xDEAD_BEEF_CAFE_BABE);
1488
1489 for _ in 0..n_reps {
1490 let o: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam / c)).collect();
1491 let s: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam * t_true)).collect();
1492 let model = ConstModel { n_e: n_bins };
1493 let obj = JointPoissonObjective {
1494 model: &model,
1495 o: &o,
1496 s: &s,
1497 c,
1498 active_mask: None,
1499 };
1500
1501 // 1D grid search over T, then local refinement via Brent-like
1502 // bisection on the gradient sign.
1503 let grid: Vec<f64> = (0..200).map(|i| 0.01 + 0.99 * (i as f64) / 199.0).collect();
1504 let mut best = (grid[0], f64::INFINITY);
1505 for &t_try in &grid {
1506 let d_try = obj
1507 .deviance_from_transmission(&vec![t_try; n_bins])
1508 .unwrap();
1509 if d_try < best.1 {
1510 best = (t_try, d_try);
1511 }
1512 }
1513 // Bisect on the gradient-sign neighbourhood.
1514 let dt = 0.01;
1515 let (mut lo, mut hi) = ((best.0 - dt).max(POISSON_EPSILON), (best.0 + dt).min(0.999));
1516 let grad_at = |t: f64| -> f64 {
1517 let tvec = vec![t; n_bins];
1518 let free_idx = [0_usize];
1519 let g = obj
1520 .deviance_gradient_analytical(&[t], &free_idx)
1521 .unwrap()
1522 .unwrap();
1523 // gradient is w.r.t. θ_0 = T (ConstModel Jacobian is 1).
1524 let _ = tvec; // silence unused
1525 g[0]
1526 };
1527 let mut glo = grad_at(lo);
1528 let mut ghi = grad_at(hi);
1529 if glo * ghi < 0.0 {
1530 for _ in 0..80 {
1531 let mid = 0.5 * (lo + hi);
1532 let gmid = grad_at(mid);
1533 if gmid * glo < 0.0 {
1534 hi = mid;
1535 ghi = gmid;
1536 } else {
1537 lo = mid;
1538 glo = gmid;
1539 }
1540 }
1541 }
1542 let t_hat = 0.5 * (lo + hi);
1543 let d_hat = obj
1544 .deviance_from_transmission(&vec![t_hat; n_bins])
1545 .unwrap();
1546 let dof = (n_bins - 1) as f64;
1547 d_per_dof_samples.push(d_hat / dof);
1548 bias_samples.push((t_hat - t_true) / t_true);
1549 }
1550
1551 let mean_dpd: f64 = d_per_dof_samples.iter().sum::<f64>() / d_per_dof_samples.len() as f64;
1552 let mean_bias: f64 = bias_samples.iter().sum::<f64>() / bias_samples.len() as f64;
1553
1554 // Under matched model, E[D]/(n-k) → 1. Tolerate [0.85, 1.15]
1555 // with n_bins=200, n_reps=30, small λ (some low-count bins).
1556 assert!(
1557 (0.85..=1.15).contains(&mean_dpd),
1558 "D/(n-k) asymptote out of band: mean={mean_dpd}"
1559 );
1560 assert!(
1561 mean_bias.abs() < 0.02,
1562 "density bias > 2%: mean={mean_bias}"
1563 );
1564 }
1565
1566 // ------------------------------------------------------------------
1567 // Edge: zero-count bin contributes 0 deviance regardless of T.
1568 // ------------------------------------------------------------------
1569 #[test]
1570 fn test_zero_counts_contribute_zero() {
1571 let model = ConstModel { n_e: 3 };
1572 let obj = JointPoissonObjective {
1573 model: &model,
1574 o: &[0.0, 10.0, 5.0],
1575 s: &[0.0, 5.0, 2.0],
1576 c: 1.5,
1577 active_mask: None,
1578 };
1579 let d_full = obj.deviance_from_transmission(&[0.6, 0.6, 0.6]).unwrap();
1580 // Drop the zero-N bin — result must be identical.
1581 let obj_reduced = JointPoissonObjective {
1582 model: &model, // same model, we just bypass the 1st bin via data
1583 o: &[10.0, 5.0],
1584 s: &[5.0, 2.0],
1585 c: 1.5,
1586 active_mask: None,
1587 };
1588 let d_reduced = obj_reduced.deviance_from_transmission(&[0.6, 0.6]).unwrap();
1589 assert!((d_full - d_reduced).abs() < 1e-12);
1590 }
1591
1592 // ------------------------------------------------------------------
1593 // FD gradient fallback agrees with analytical form.
1594 // ------------------------------------------------------------------
1595 #[test]
1596 fn test_fd_gradient_matches_analytical() {
1597 let e: Vec<f64> = (0..15).map(|i| 0.2 + 0.1 * i as f64).collect();
1598 let theta = [0.9_f64, 0.05_f64];
1599 let c = 1.5;
1600 let lam = 300.0;
1601 let model = LinearModel { e: &e };
1602 let t_true = model.evaluate(&theta).unwrap();
1603 let o: Vec<f64> = t_true.iter().map(|_| lam / c).collect();
1604 let s: Vec<f64> = t_true.iter().map(|&ti| lam * ti).collect();
1605 let obj = JointPoissonObjective {
1606 model: &model,
1607 o: &o,
1608 s: &s,
1609 c,
1610 active_mask: None,
1611 };
1612 let mut ps = ParameterSet::new(vec![
1613 FitParameter::non_negative("theta_0", 0.85),
1614 FitParameter::non_negative("theta_1", 0.06),
1615 ]);
1616 let g_fd = obj.deviance_gradient_fd(&mut ps, 1e-6).unwrap();
1617 let g_analytical = obj
1618 .deviance_gradient_analytical(&ps.all_values(), &ps.free_indices())
1619 .unwrap()
1620 .unwrap();
1621 for (f, a) in g_fd.iter().zip(g_analytical.iter()) {
1622 let rel = ((f - a) / a.abs().max(1e-6)).abs();
1623 assert!(rel < 5e-3, "fd={f} analytical={a} rel={rel}");
1624 }
1625 }
1626
1627 // ------------------------------------------------------------------
1628 // Fisher matrix is symmetric positive semi-definite at the fit.
1629 // ------------------------------------------------------------------
1630 #[test]
1631 fn test_fisher_matrix_symmetry_psd() {
1632 let e: Vec<f64> = (0..10).map(|i| 0.3 + 0.1 * i as f64).collect();
1633 let theta = [0.9_f64, 0.05_f64];
1634 let c = 2.0;
1635 let lam = 400.0;
1636 let model = LinearModel { e: &e };
1637 let t_true = model.evaluate(&theta).unwrap();
1638 let o: Vec<f64> = t_true.iter().map(|_| lam / c).collect();
1639 let s: Vec<f64> = t_true.iter().map(|&ti| lam * ti).collect();
1640 let obj = JointPoissonObjective {
1641 model: &model,
1642 o: &o,
1643 s: &s,
1644 c,
1645 active_mask: None,
1646 };
1647 let info = obj
1648 .fisher_information(&theta, &[0, 1])
1649 .unwrap()
1650 .expect("LinearModel provides analytical jacobian");
1651 // Symmetry.
1652 let i01 = info.get(0, 1);
1653 let i10 = info.get(1, 0);
1654 assert!((i01 - i10).abs() < 1e-10);
1655 // PSD: diagonal entries > 0 (model is identifiable).
1656 assert!(info.get(0, 0) > 0.0);
1657 assert!(info.get(1, 1) > 0.0);
1658 // Determinant > 0 (rank-2 identifiable).
1659 let det = info.get(0, 0) * info.get(1, 1) - i01 * i10;
1660 assert!(det > 0.0, "Fisher matrix determinant = {det}");
1661 }
1662
1663 // ==================================================================
1664 // joint_poisson_fit — end-to-end integration tests
1665 // ==================================================================
1666
1667 /// A wrapped transmission model: T_out = A_n · T_inner + B_A + B_B/√E + B_C·√E.
1668 /// Models the full counts-path background structure of memo 35 §P2.2.
1669 struct BackgroundedTransmission<'a> {
1670 inner: &'a dyn FitModel,
1671 energies: &'a [f64],
1672 n_idx: usize,
1673 a_idx: usize,
1674 b_a_idx: usize,
1675 b_b_idx: usize,
1676 b_c_idx: usize,
1677 n_params: usize,
1678 }
1679
1680 impl<'a> FitModel for BackgroundedTransmission<'a> {
1681 fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1682 // Pass the "density" parameter to the inner model as its param 0.
1683 let t_inner = self.inner.evaluate(&[params[self.n_idx]])?;
1684 let a_n = params[self.a_idx];
1685 let b_a = params[self.b_a_idx];
1686 let b_b = params[self.b_b_idx];
1687 let b_c = params[self.b_c_idx];
1688 Ok(t_inner
1689 .iter()
1690 .zip(self.energies.iter())
1691 .map(|(&t, &e)| {
1692 let inv_sqrt_e = if e > 0.0 { 1.0 / e.sqrt() } else { 0.0 };
1693 let sqrt_e = if e > 0.0 { e.sqrt() } else { 0.0 };
1694 a_n * t + b_a + b_b * inv_sqrt_e + b_c * sqrt_e
1695 })
1696 .collect())
1697 }
1698 // No analytical jacobian — forces the fitter onto FD fallback, which
1699 // is the stress test (memo 35 §P2.1 notes FD + over-parameterization
1700 // as the stall trigger).
1701 }
1702
1703 /// Exponential-in-E model: T_inner = exp(−n · σ(E)), σ(E) = 1.
1704 /// Effectively a single-parameter constant transmission when σ=1 flat.
1705 /// Uses an energy-dependent "cross section" so Jacobian is identifiable.
1706 struct ExpDecayModel<'a> {
1707 sigma: &'a [f64],
1708 }
1709 impl<'a> FitModel for ExpDecayModel<'a> {
1710 fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1711 let n = params[0];
1712 Ok(self
1713 .sigma
1714 .iter()
1715 .map(|&s| (-n * s).exp().max(POISSON_EPSILON))
1716 .collect())
1717 }
1718 fn analytical_jacobian(
1719 &self,
1720 _params: &[f64],
1721 free_param_indices: &[usize],
1722 y_current: &[f64],
1723 ) -> Option<FlatMatrix> {
1724 // ∂T/∂n = -σ · T
1725 let n_e = y_current.len();
1726 let n_free = free_param_indices.len();
1727 let mut jac = FlatMatrix::zeros(n_e, n_free);
1728 for (i, &y_i) in y_current.iter().enumerate() {
1729 for (j, &pi) in free_param_indices.iter().enumerate() {
1730 *jac.get_mut(i, j) = if pi == 0 { -self.sigma[i] * y_i } else { 0.0 };
1731 }
1732 }
1733 Some(jac)
1734 }
1735 }
1736
1737 /// Deterministic Poisson generator (Knuth for small λ, Gaussian for
1738 /// large). Shared across the stochastic-asymptote and joint-Poisson
1739 /// fit tests in this module.
1740 struct Xorshift(u64);
1741 impl Xorshift {
1742 fn next_u64(&mut self) -> u64 {
1743 let mut x = self.0;
1744 x ^= x << 13;
1745 x ^= x >> 7;
1746 x ^= x << 17;
1747 self.0 = x;
1748 x
1749 }
1750 fn uniform(&mut self) -> f64 {
1751 (self.next_u64() as f64) / (u64::MAX as f64)
1752 }
1753 fn poisson(&mut self, lambda: f64) -> f64 {
1754 if lambda <= 0.0 {
1755 return 0.0;
1756 }
1757 if lambda > 30.0 {
1758 let u1 = self.uniform().max(1e-12);
1759 let u2 = self.uniform();
1760 let z = (-2.0 * u1.ln()).sqrt() * (2.0 * std::f64::consts::PI * u2).cos();
1761 return (lambda + z * lambda.sqrt()).round().max(0.0);
1762 }
1763 let l = (-lambda).exp();
1764 let mut k: f64 = 0.0;
1765 let mut p: f64 = 1.0;
1766 loop {
1767 k += 1.0;
1768 let u = self.uniform();
1769 p *= u;
1770 if p <= l {
1771 return k - 1.0;
1772 }
1773 if k > 1000.0 {
1774 return k - 1.0;
1775 }
1776 }
1777 }
1778 }
1779
1780 // ------------------------------------------------------------------
1781 // Matched-model single-parameter recovery at c = 5.98.
1782 // This is the EG1 "proposed" cell in miniature — verify |bias| < 1%
1783 // and D / (n − k) ∈ [0.85, 1.15] without needing the polish.
1784 // ------------------------------------------------------------------
1785 #[test]
1786 fn test_joint_poisson_fit_matched_model_single_param() {
1787 // Energies 1..10, flat cross section σ = 1. Truth n = 0.3.
1788 let n_bins = 200;
1789 let sigma = vec![1.0_f64; n_bins];
1790 let model = ExpDecayModel { sigma: &sigma };
1791 let n_true = 0.3_f64;
1792 let c = 5.98;
1793 let lam = 3000.0; // OB target ~500 counts/bin
1794 let t_true = model.evaluate(&[n_true]).unwrap();
1795
1796 let mut rng = Xorshift(0x1234_5678_9ABC_DEF0);
1797 let o: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam / c)).collect();
1798 let s: Vec<f64> = (0..n_bins).map(|i| rng.poisson(lam * t_true[i])).collect();
1799
1800 let obj = JointPoissonObjective {
1801 model: &model,
1802 o: &o,
1803 s: &s,
1804 c,
1805 active_mask: None,
1806 };
1807 let mut params = ParameterSet::new(vec![FitParameter::non_negative("n", 0.1)]);
1808 let cfg = JointPoissonFitConfig {
1809 enable_polish: true,
1810 ..Default::default()
1811 };
1812 let result = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
1813
1814 let n_fit = result.params[0];
1815 let rel_bias = (n_fit - n_true) / n_true;
1816 assert!(
1817 rel_bias.abs() < 0.01,
1818 "density bias {rel_bias} exceeds 1% (n_fit={n_fit} n_true={n_true})"
1819 );
1820 assert!(
1821 (0.85..=1.15).contains(&result.deviance_per_dof),
1822 "D/(n-k) out of band: {}",
1823 result.deviance_per_dof
1824 );
1825 }
1826
1827 // ------------------------------------------------------------------
1828 // Polish-never-worsens invariant on a backgrounded fit. Memo 35 §P2.1
1829 // claims NM polish reduces D materially when stage-1 stalls. At the
1830 // unit-test scale we verify the testable invariant: enabling polish
1831 // never produces a larger final D than disabling it on the same data.
1832 //
1833 // Note: on this over-parameterized (5-free-param) synthetic with only
1834 // 150 bins, the deviance surface has multiple near-equal minima —
1835 // exactly the identifiability ambiguity §P2.2 targets. Density
1836 // recovery under over-parameterization is therefore *not* a unit-test
1837 // contract here; it is tested end-to-end with the single-parameter
1838 // matched-model test above.
1839 // ------------------------------------------------------------------
1840 #[test]
1841 fn test_joint_poisson_fit_polish_does_not_worsen_deviance() {
1842 let n_bins = 150;
1843 let energies: Vec<f64> = (0..n_bins).map(|i| 1.0 + 0.5 * i as f64).collect();
1844 let sigma: Vec<f64> = energies.iter().map(|&e| 1.0 / e).collect();
1845 let inner = ExpDecayModel { sigma: &sigma };
1846
1847 // Truth: n = 0.3, A_n = 0.9, no additive bg.
1848 let n_true = 0.3_f64;
1849 let a_n_true = 0.9_f64;
1850 let t_inner_true = inner.evaluate(&[n_true]).unwrap();
1851 let t_true: Vec<f64> = t_inner_true.iter().map(|&t| a_n_true * t).collect();
1852
1853 let c = 5.98_f64;
1854 let lam = 5000.0_f64;
1855 let mut rng = Xorshift(0xF00D_FACE_DEAD_BEEF);
1856 let o: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam / c)).collect();
1857 let s: Vec<f64> = (0..n_bins).map(|i| rng.poisson(lam * t_true[i])).collect();
1858
1859 let bg_model = BackgroundedTransmission {
1860 inner: &inner,
1861 energies: &energies,
1862 n_idx: 0,
1863 a_idx: 1,
1864 b_a_idx: 2,
1865 b_b_idx: 3,
1866 b_c_idx: 4,
1867 n_params: 5,
1868 };
1869 let _ = bg_model.n_params; // silence dead-code warning
1870
1871 let obj = JointPoissonObjective {
1872 model: &bg_model,
1873 o: &o,
1874 s: &s,
1875 c,
1876 active_mask: None,
1877 };
1878
1879 // x0 analogous to EG2-S1 regime: n near truth, A_n = 1, all
1880 // additive bg at 0, bg bounds tight to curb degeneracy.
1881 let mk_params = || {
1882 ParameterSet::new(vec![
1883 FitParameter::non_negative("n", 0.25),
1884 FitParameter::non_negative("A_n", 1.0),
1885 FitParameter {
1886 name: "B_A".into(),
1887 value: 0.0,
1888 lower: -0.05,
1889 upper: 0.05,
1890 fixed: false,
1891 },
1892 FitParameter {
1893 name: "B_B".into(),
1894 value: 0.0,
1895 lower: -0.05,
1896 upper: 0.05,
1897 fixed: false,
1898 },
1899 FitParameter {
1900 name: "B_C".into(),
1901 value: 0.0,
1902 lower: -0.05,
1903 upper: 0.05,
1904 fixed: false,
1905 },
1906 ])
1907 };
1908
1909 let mut params_no_polish = mk_params();
1910 let cfg_no_polish = JointPoissonFitConfig {
1911 enable_polish: false,
1912 ..Default::default()
1913 };
1914 let r_no_polish = joint_poisson_fit(&obj, &mut params_no_polish, &cfg_no_polish).unwrap();
1915
1916 let mut params_polish = mk_params();
1917 let cfg_polish = JointPoissonFitConfig {
1918 enable_polish: true,
1919 ..Default::default()
1920 };
1921 let r_polish = joint_poisson_fit(&obj, &mut params_polish, &cfg_polish).unwrap();
1922
1923 // Invariant: enabling polish must not increase final D.
1924 assert!(
1925 r_polish.deviance <= r_no_polish.deviance + 1e-6,
1926 "polish worsened D: D_polish={} D_no_polish={}",
1927 r_polish.deviance,
1928 r_no_polish.deviance
1929 );
1930
1931 // When polish_improved flag is set, polish D must be strictly
1932 // better than stage-1 D (consistency check on the flag semantics).
1933 if r_polish.polish_improved {
1934 assert!(
1935 r_polish.deviance < r_no_polish.deviance,
1936 "polish_improved=true but D_polish={} >= D_no_polish={}",
1937 r_polish.deviance,
1938 r_no_polish.deviance
1939 );
1940 }
1941
1942 // The fit should return a physically sensible density (positive,
1943 // finite, within an order of magnitude of truth — not a strict
1944 // recovery test, just a sanity check).
1945 let n_fit = r_polish.params[0];
1946 assert!(n_fit.is_finite() && n_fit > 0.0);
1947 assert!(
1948 n_fit > 0.1 && n_fit < 0.8,
1949 "density grossly off: n_fit={n_fit} (truth={n_true})"
1950 );
1951 }
1952
1953 // ------------------------------------------------------------------
1954 // Fit result carries gn_converged and polish_converged separately
1955 // (memo 35 §P2.3 — acceptance from deviance value, not one flag).
1956 // ------------------------------------------------------------------
1957 #[test]
1958 fn test_joint_poisson_fit_exposes_separate_converged_flags() {
1959 let n_bins = 50;
1960 let sigma = vec![0.5_f64; n_bins];
1961 let model = ExpDecayModel { sigma: &sigma };
1962 let n_true = 0.2;
1963 let c = 2.0;
1964 let lam = 500.0;
1965 let t_true = model.evaluate(&[n_true]).unwrap();
1966 let mut rng = Xorshift(0xABAD_CAFE_BABE_F00D);
1967 let o: Vec<f64> = (0..n_bins).map(|_| rng.poisson(lam / c)).collect();
1968 let s: Vec<f64> = (0..n_bins).map(|i| rng.poisson(lam * t_true[i])).collect();
1969
1970 let obj = JointPoissonObjective {
1971 model: &model,
1972 o: &o,
1973 s: &s,
1974 c,
1975 active_mask: None,
1976 };
1977 let mut params = ParameterSet::new(vec![FitParameter::non_negative("n", 0.1)]);
1978 let cfg = JointPoissonFitConfig {
1979 enable_polish: true,
1980 ..Default::default()
1981 };
1982 let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
1983
1984 // Both flags exist; at least one should be true on this easy case.
1985 assert!(r.gn_converged || r.polish_converged);
1986 assert!(r.n_data == n_bins);
1987 assert!(r.n_free == 1);
1988 assert!(r.deviance > 0.0);
1989 assert!(r.deviance_per_dof.is_finite());
1990 // Uncertainty present (compute_covariance default true).
1991 assert!(r.uncertainties.is_some());
1992 let u = r.uncertainties.as_ref().unwrap();
1993 assert_eq!(u.len(), 1);
1994 assert!(u[0].is_finite() && u[0] > 0.0);
1995 }
1996
1997 // ------------------------------------------------------------------
1998 // Reported uncertainty matches the analytical Cramér-Rao bound
1999 // I^{-1} (NOT (2I)^{-1} — the Hessian-of-D inverse, which would
2000 // under-report σ by √2). Caught in code review of memo-35 §P1
2001 // implementation; see `joint_poisson_fit` covariance-extraction
2002 // doc-comment for the rescaling rationale.
2003 // ------------------------------------------------------------------
2004 #[test]
2005 fn test_uncertainty_matches_analytical_fisher_inverse() {
2006 // Construct a single-parameter constant-T model on noise-free
2007 // expected counts: O_i = λ/c, S_i = λ·T (per memo 35 §4.1).
2008 // With ConstModel (J_i = ∂T/∂θ = 1), the analytical Fisher is
2009 // I(T) = Σ_i (O_i + S_i)·c / (T·(1+cT)²)
2010 // = N · λ · (1+cT)/c · c / (T·(1+cT)²)
2011 // = N · λ / (T · (1+cT))
2012 // and σ_T = √(I^{-1}) = √( T·(1+cT) / (N·λ) ).
2013 let n_bins = 200;
2014 let t_true = 0.5_f64;
2015 let c = 2.0_f64;
2016 let lam = 100.0_f64;
2017 let o: Vec<f64> = vec![lam / c; n_bins];
2018 let s: Vec<f64> = vec![lam * t_true; n_bins];
2019 let model = ConstModel { n_e: n_bins };
2020 let obj = JointPoissonObjective {
2021 model: &model,
2022 o: &o,
2023 s: &s,
2024 c,
2025 active_mask: None,
2026 };
2027 let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", t_true)]);
2028 let cfg = JointPoissonFitConfig {
2029 // Disable polish for a clean Newton-only fit (avoids NM-tail
2030 // perturbations of the final θ that would shift σ slightly).
2031 enable_polish: false,
2032 ..Default::default()
2033 };
2034 let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
2035 let sigma_reported = r.uncertainties.as_ref().expect("σ available")[0];
2036
2037 // Analytical Cramér-Rao σ.
2038 let sigma_analytical = (t_true * (1.0 + c * t_true) / (n_bins as f64 * lam)).sqrt();
2039
2040 // The pre-fix (uncompensated) value would be σ_analytical / √2 —
2041 // tighten the tolerance below √2 so the regression is caught.
2042 let rel_err = (sigma_reported - sigma_analytical).abs() / sigma_analytical;
2043 assert!(
2044 rel_err < 0.05,
2045 "reported σ = {sigma_reported} vs analytical I^{{-1}}^(1/2) = \
2046 {sigma_analytical} (rel_err = {rel_err}); pre-fix code reported \
2047 σ_analytical / √2 ≈ {} which would give rel_err ≈ 0.293",
2048 sigma_analytical / 2.0_f64.sqrt(),
2049 );
2050 }
2051
2052 // ------------------------------------------------------------------
2053 // Active-bin mask (SAMMY EMIN/EMAX-equivalent fit-energy-range, #514).
2054 // ------------------------------------------------------------------
2055
2056 /// `deviance_from_transmission` with `active_mask` set must equal
2057 /// the same call computed only over the `true` bins (subset
2058 /// equivalence) — the masking is correct iff dropping out-of-mask
2059 /// bins from `o`, `s`, `t` produces the same value.
2060 #[test]
2061 fn test_jp_active_mask_subset_equivalence() {
2062 // 5-bin objective with an arbitrary mask — bins 1 and 3 active.
2063 let o_full = [10.0, 20.0, 5.0, 15.0, 25.0];
2064 let s_full = [4.0, 8.0, 1.0, 6.0, 12.0];
2065 let t_full = [0.4, 0.5, 0.7, 0.6, 0.45];
2066 let mask = [false, true, false, true, false];
2067 let c = 1.5;
2068 let model_full = ConstModel { n_e: 5 };
2069 let obj_full = JointPoissonObjective {
2070 model: &model_full,
2071 o: &o_full,
2072 s: &s_full,
2073 c,
2074 active_mask: Some(&mask),
2075 };
2076 let d_masked = obj_full.deviance_from_transmission(&t_full).unwrap();
2077
2078 // Compare against an objective built directly on the active subset.
2079 let o_sub = [o_full[1], o_full[3]];
2080 let s_sub = [s_full[1], s_full[3]];
2081 let t_sub = [t_full[1], t_full[3]];
2082 let model_sub = ConstModel { n_e: 2 };
2083 let obj_sub = JointPoissonObjective {
2084 model: &model_sub,
2085 o: &o_sub,
2086 s: &s_sub,
2087 c,
2088 active_mask: None,
2089 };
2090 let d_subset = obj_sub.deviance_from_transmission(&t_sub).unwrap();
2091
2092 assert!(
2093 (d_masked - d_subset).abs() < 1e-12,
2094 "masked deviance {d_masked} != subset deviance {d_subset}"
2095 );
2096
2097 // Active-bin count should be 2, not 5.
2098 assert_eq!(obj_full.n_active(), 2);
2099 assert_eq!(obj_full.n_data(), 5);
2100 }
2101
2102 /// Out-of-mask gradient contributions must drop to zero — verified
2103 /// by comparing against an unmasked subset gradient.
2104 #[test]
2105 fn test_jp_active_mask_gradient_subset_equivalence() {
2106 let e_full: Vec<f64> = (0..6).map(|i| 0.1 + 0.1 * i as f64).collect();
2107 let theta_true = [0.95_f64, 0.05_f64];
2108 let c = 2.0;
2109 let lam = 100.0;
2110 let model_full = LinearModel { e: &e_full };
2111 let t_full = model_full.evaluate(&theta_true).unwrap();
2112 let o_full: Vec<f64> = vec![lam / c; e_full.len()];
2113 let s_full: Vec<f64> = t_full.iter().map(|&ti| lam * ti).collect();
2114
2115 // Mask = bins 2..5 active.
2116 let mask = vec![false, false, true, true, true, false];
2117 let obj_full = JointPoissonObjective {
2118 model: &model_full,
2119 o: &o_full,
2120 s: &s_full,
2121 c,
2122 active_mask: Some(&mask),
2123 };
2124
2125 let params_full = ParameterSet::new(vec![
2126 FitParameter::non_negative("a", theta_true[0]),
2127 FitParameter::non_negative("b", theta_true[1]),
2128 ]);
2129 let free_idx = params_full.free_indices();
2130 let theta_eval = [0.9_f64, 0.07_f64];
2131 let grad_masked = obj_full
2132 .deviance_gradient_analytical(&theta_eval, &free_idx)
2133 .unwrap()
2134 .expect("analytical gradient");
2135
2136 // Subset reference: only bins 2..5.
2137 let e_sub = e_full[2..5].to_vec();
2138 let o_sub = o_full[2..5].to_vec();
2139 let s_sub = s_full[2..5].to_vec();
2140 let model_sub = LinearModel { e: &e_sub };
2141 let obj_sub = JointPoissonObjective {
2142 model: &model_sub,
2143 o: &o_sub,
2144 s: &s_sub,
2145 c,
2146 active_mask: None,
2147 };
2148 let grad_sub = obj_sub
2149 .deviance_gradient_analytical(&theta_eval, &free_idx)
2150 .unwrap()
2151 .expect("analytical gradient");
2152
2153 for (i, (&gm, &gs)) in grad_masked.iter().zip(grad_sub.iter()).enumerate() {
2154 assert!(
2155 (gm - gs).abs() < 1e-9,
2156 "grad component {i}: masked={gm} subset={gs}"
2157 );
2158 }
2159 }
2160
2161 /// `joint_poisson_fit` must reject an underdetermined (n_active <
2162 /// n_free) configuration with a non-converged result and NaN
2163 /// deviance / per-dof, mirroring the LM solver. An all-`false`
2164 /// active mask is the extreme case (`n_active == 0 < n_free`);
2165 /// the prior `.max(1)` divisor produced a deceptive
2166 /// finite-looking deviance-per-dof for empty / too-narrow masks.
2167 /// Regression for Round-2 review fix #2 (#514).
2168 #[test]
2169 fn test_joint_poisson_rejects_zero_active_mask() {
2170 let n_bins = 10;
2171 let o: Vec<f64> = vec![50.0; n_bins];
2172 let s: Vec<f64> = vec![25.0; n_bins];
2173 let mask = vec![false; n_bins]; // n_active = 0
2174 let model = ConstModel { n_e: n_bins };
2175 let obj = JointPoissonObjective {
2176 model: &model,
2177 o: &o,
2178 s: &s,
2179 c: 1.0,
2180 active_mask: Some(&mask),
2181 };
2182 let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
2183 let cfg = JointPoissonFitConfig::default();
2184 let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
2185
2186 assert!(
2187 !r.gn_converged && !r.polish_converged,
2188 "underdetermined fit must report non-converged"
2189 );
2190 assert!(
2191 r.deviance.is_nan(),
2192 "underdetermined deviance must be NaN, got {}",
2193 r.deviance
2194 );
2195 assert!(
2196 r.deviance_per_dof.is_nan(),
2197 "underdetermined deviance-per-dof must be NaN, got {}",
2198 r.deviance_per_dof
2199 );
2200 assert_eq!(r.n_data, n_bins);
2201 assert_eq!(r.n_active, 0);
2202 assert_eq!(r.n_free, 1);
2203 assert!(r.covariance.is_none());
2204 assert!(r.uncertainties.is_none());
2205 }
2206
2207 /// Zero active bins with **all parameters fixed** (`n_free == 0`)
2208 /// must still return non-converged. Without the
2209 /// `n_active == 0` early-return, the underdetermined check
2210 /// `n_active < n_free` is `0 < 0` → false, so the function would
2211 /// fall through to the main loop, compute `deviance = 0` from the
2212 /// empty sum, and `dof = 0` → `deviance_per_dof = NaN` — but
2213 /// `gn_converged` could still be `true`, masquerading as a
2214 /// successful fit on no data. Regression for #517 R3 P1 (#514).
2215 #[test]
2216 fn test_joint_poisson_rejects_zero_active_with_no_free_params() {
2217 let n_bins = 5;
2218 let o: Vec<f64> = vec![10.0; n_bins];
2219 let s: Vec<f64> = vec![5.0; n_bins];
2220 let mask = vec![false; n_bins];
2221 let model = ConstModel { n_e: n_bins };
2222 let obj = JointPoissonObjective {
2223 model: &model,
2224 o: &o,
2225 s: &s,
2226 c: 1.0,
2227 active_mask: Some(&mask),
2228 };
2229 let mut params = ParameterSet::new(vec![FitParameter::fixed("T", 0.5)]);
2230 let r = joint_poisson_fit(&obj, &mut params, &JointPoissonFitConfig::default()).unwrap();
2231 assert!(!r.gn_converged);
2232 assert!(!r.polish_converged);
2233 assert!(r.deviance.is_nan());
2234 assert!(r.deviance_per_dof.is_nan());
2235 assert_eq!(r.n_active, 0);
2236 assert_eq!(r.n_free, 0);
2237 }
2238
2239 /// `joint_poisson_fit` validates active-mask length up-front and
2240 /// returns `LengthMismatch` rather than relying on a debug-assert
2241 /// deep in the deviance routines (which silently passes through in
2242 /// release builds, then panics on out-of-bounds index reads).
2243 /// Regression for Round-2 review fix #5 (#514).
2244 #[test]
2245 fn test_joint_poisson_rejects_active_mask_length_mismatch() {
2246 let n_bins = 5;
2247 let o: Vec<f64> = vec![10.0; n_bins];
2248 let s: Vec<f64> = vec![5.0; n_bins];
2249 let mask_wrong = vec![true, true, true]; // wrong length
2250 let model = ConstModel { n_e: n_bins };
2251 let obj = JointPoissonObjective {
2252 model: &model,
2253 o: &o,
2254 s: &s,
2255 c: 1.0,
2256 active_mask: Some(&mask_wrong),
2257 };
2258 let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
2259 let cfg = JointPoissonFitConfig::default();
2260 let err = joint_poisson_fit(&obj, &mut params, &cfg).unwrap_err();
2261 assert!(
2262 matches!(
2263 err,
2264 FittingError::LengthMismatch {
2265 field: "active_mask",
2266 ..
2267 }
2268 ),
2269 "expected LengthMismatch on active_mask; got {err:?}"
2270 );
2271 }
2272
2273 // ==================================================================
2274 // Release-mode input validation at joint_poisson_fit.
2275 //
2276 // The inner `binomial_deviance_term` and `deviance_from_transmission`
2277 // protect themselves with `debug_assert!` only. Release builds skip
2278 // those, so a length mismatch in `o` vs `s` silently truncates via
2279 // `.zip()` and a non-positive `c` produces finite garbage that the
2280 // optimizer happily minimises. Validate at the public entry point.
2281 // ==================================================================
2282
2283 /// `joint_poisson_fit` rejects an `o`/`s` length mismatch with a
2284 /// `LengthMismatch` error rather than silently truncating via `.zip()`
2285 /// and minimising bogus deviance on a sub-range of bins.
2286 #[test]
2287 fn test_joint_poisson_rejects_o_s_length_mismatch() {
2288 let n_bins = 5;
2289 let o: Vec<f64> = vec![10.0; n_bins];
2290 // Deliberate mismatch: `s` has one fewer bin than `o`.
2291 let s: Vec<f64> = vec![5.0; n_bins - 1];
2292 let model = ConstModel { n_e: n_bins };
2293 let obj = JointPoissonObjective {
2294 model: &model,
2295 o: &o,
2296 s: &s,
2297 c: 1.0,
2298 active_mask: None,
2299 };
2300 let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
2301 let err =
2302 joint_poisson_fit(&obj, &mut params, &JointPoissonFitConfig::default()).unwrap_err();
2303 assert!(
2304 matches!(
2305 err,
2306 FittingError::LengthMismatch {
2307 field: "sample_counts",
2308 ..
2309 }
2310 ),
2311 "expected LengthMismatch on sample_counts; got {err:?}"
2312 );
2313 }
2314
2315 /// `joint_poisson_fit` rejects a non-positive proton-charge ratio `c`
2316 /// with `InvalidConfig` rather than falling through to the inner
2317 /// `debug_assert!` (which is a no-op in release builds and lets the
2318 /// optimizer minimise a garbage deviance landscape).
2319 #[test]
2320 fn test_joint_poisson_rejects_non_positive_c() {
2321 let n_bins = 5;
2322 let o: Vec<f64> = vec![10.0; n_bins];
2323 let s: Vec<f64> = vec![5.0; n_bins];
2324 let model = ConstModel { n_e: n_bins };
2325 let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
2326 // c = 0 is the textbook degenerate case (no sample counts).
2327 let obj_zero = JointPoissonObjective {
2328 model: &model,
2329 o: &o,
2330 s: &s,
2331 c: 0.0,
2332 active_mask: None,
2333 };
2334 let err = joint_poisson_fit(&obj_zero, &mut params, &JointPoissonFitConfig::default())
2335 .unwrap_err();
2336 assert!(
2337 matches!(err, FittingError::InvalidConfig(_)),
2338 "expected InvalidConfig on c=0; got {err:?}"
2339 );
2340
2341 // Negative c is unphysical.
2342 let mut params2 = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
2343 let obj_neg = JointPoissonObjective {
2344 model: &model,
2345 o: &o,
2346 s: &s,
2347 c: -1.5,
2348 active_mask: None,
2349 };
2350 let err = joint_poisson_fit(&obj_neg, &mut params2, &JointPoissonFitConfig::default())
2351 .unwrap_err();
2352 assert!(
2353 matches!(err, FittingError::InvalidConfig(_)),
2354 "expected InvalidConfig on c<0; got {err:?}"
2355 );
2356
2357 // NaN c — caught by the same finiteness check.
2358 let mut params3 = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
2359 let obj_nan = JointPoissonObjective {
2360 model: &model,
2361 o: &o,
2362 s: &s,
2363 c: f64::NAN,
2364 active_mask: None,
2365 };
2366 let err = joint_poisson_fit(&obj_nan, &mut params3, &JointPoissonFitConfig::default())
2367 .unwrap_err();
2368 assert!(
2369 matches!(err, FittingError::InvalidConfig(_)),
2370 "expected InvalidConfig on c=NaN; got {err:?}"
2371 );
2372 }
2373
2374 // ==================================================================
2375 // `f64::max(NaN, ε) == ε` swallows active NaN T.
2376 //
2377 // Rust stdlib's `f64::max` returns the non-NaN argument when one is
2378 // NaN, so `t.max(POISSON_EPSILON)` silently turns a NaN transmission
2379 // into ε. The deviance term then evaluates to a finite (large)
2380 // number which passes the trial-step's `v.is_finite()` guard, so the
2381 // optimizer accepts steps into regions where the model is broken.
2382 //
2383 // `binomial_deviance_term` returns NaN when T is non-finite (so the
2384 // deviance sum becomes NaN and the trial guard rejects the step),
2385 // and `deviance_weight` / `deviance_curvature` return 0 (so the
2386 // gradient / Fisher accumulators are not poisoned by the bad bin).
2387 // ==================================================================
2388
2389 /// `binomial_deviance_term` returns NaN when `t` is non-finite — so
2390 /// the per-bin sum poisons the deviance and the trial-step guard
2391 /// (`Ok(v) if v.is_finite()`) rejects the step instead of silently
2392 /// accepting a bogus-but-finite value.
2393 #[test]
2394 fn test_binomial_deviance_term_nan_t_returns_nan() {
2395 // Pre-fix: `t.max(POISSON_EPSILON)` swallows NaN and returns a
2396 // finite (but meaningless) deviance.
2397 let d_nan_t = binomial_deviance_term(50.0, 10.0, f64::NAN, 2.0);
2398 assert!(
2399 d_nan_t.is_nan(),
2400 "non-finite T must produce NaN deviance, not a finite shim; got {d_nan_t}"
2401 );
2402
2403 // +inf / -inf likewise — they are not physical transmission values.
2404 let d_inf_t = binomial_deviance_term(50.0, 10.0, f64::INFINITY, 2.0);
2405 assert!(
2406 d_inf_t.is_nan(),
2407 "+inf T must produce NaN deviance; got {d_inf_t}"
2408 );
2409 let d_neg_inf_t = binomial_deviance_term(50.0, 10.0, f64::NEG_INFINITY, 2.0);
2410 assert!(
2411 d_neg_inf_t.is_nan(),
2412 "-inf T must produce NaN deviance; got {d_neg_inf_t}"
2413 );
2414 }
2415
2416 /// `deviance_weight` returns 0 for non-finite `t` so the gradient
2417 /// accumulator is not poisoned — bad bins drop out instead of
2418 /// becoming silent NaN contributions weighted by the Jacobian.
2419 #[test]
2420 fn test_deviance_weight_nan_t_returns_zero() {
2421 let w = deviance_weight(50.0, 10.0, f64::NAN, 2.0);
2422 assert_eq!(w, 0.0, "non-finite T must give zero weight; got {w}");
2423 }
2424
2425 /// `deviance_curvature` returns 0 for non-finite `t` so the Fisher
2426 /// info accumulator is not poisoned.
2427 #[test]
2428 fn test_deviance_curvature_nan_t_returns_zero() {
2429 let h = deviance_curvature(50.0, 10.0, f64::NAN, 2.0);
2430 assert_eq!(h, 0.0, "non-finite T must give zero curvature; got {h}");
2431 }
2432
2433 /// End-to-end: a model that returns NaN at some active bin makes the
2434 /// deviance non-finite, the trial-step guard rejects it (rather than
2435 /// accepting a bogus finite step), and the fit either bails out
2436 /// non-converged or recovers without committing the bad step. Prior
2437 /// to the M14 fix the optimizer could silently accept the NaN step.
2438 #[test]
2439 fn test_joint_poisson_fit_rejects_nan_transmission() {
2440 // Model that returns NaN at θ < 0.1 and a constant 0.5 otherwise.
2441 struct NanAtSmallTheta;
2442 impl FitModel for NanAtSmallTheta {
2443 fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
2444 let t = if params[0] < 0.1 { f64::NAN } else { 0.5 };
2445 Ok(vec![t; 4])
2446 }
2447 fn analytical_jacobian(
2448 &self,
2449 _params: &[f64],
2450 free_param_indices: &[usize],
2451 y_current: &[f64],
2452 ) -> Option<FlatMatrix> {
2453 let n_e = y_current.len();
2454 let n_free = free_param_indices.len();
2455 let mut jac = FlatMatrix::zeros(n_e, n_free);
2456 for i in 0..n_e {
2457 for (j, &pi) in free_param_indices.iter().enumerate() {
2458 *jac.get_mut(i, j) = if pi == 0 { 1.0 } else { 0.0 };
2459 }
2460 }
2461 Some(jac)
2462 }
2463 }
2464
2465 let model = NanAtSmallTheta;
2466 let n = 4;
2467 let o = vec![10.0; n];
2468 let s = vec![5.0; n];
2469 let obj = JointPoissonObjective {
2470 model: &model,
2471 o: &o,
2472 s: &s,
2473 c: 1.0,
2474 active_mask: None,
2475 };
2476 // Initial point lands in the NaN region.
2477 let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.05)]);
2478 let cfg = JointPoissonFitConfig::default();
2479 let result = joint_poisson_fit(&obj, &mut params, &cfg);
2480 match result {
2481 Ok(r) => {
2482 // The optimizer must NOT report a finite deviance from a
2483 // NaN-T initial point — pre-fix it would do so by silently
2484 // converting NaN to POISSON_EPSILON. After the fix the
2485 // deviance is NaN (initial eval propagates), or the fit
2486 // never accepts a NaN step, or it ascends out of the NaN
2487 // region and lands at the finite plateau (params[0] >= 0.1).
2488 if r.params[0] < 0.1 {
2489 assert!(
2490 r.deviance.is_nan() && !r.gn_converged,
2491 "stayed in NaN region but reported finite deviance: {r:?}"
2492 );
2493 }
2494 }
2495 Err(_) => {
2496 // Acceptable: hard error from the initial evaluation.
2497 }
2498 }
2499 }
2500
2501 /// All-fixed parameters + NaN transmission must NOT be reported as
2502 /// `gn_converged = true`.
2503 ///
2504 /// The `n_free == 0` shortcut in `damped_fisher_stage` previously set
2505 /// `converged = true` unconditionally, so a fit with every parameter
2506 /// fixed and a model that returns NaN at active bins would return
2507 /// `deviance = NaN` together with `gn_converged = true`. Downstream
2508 /// pipeline code (`pipeline.rs`'s `gn_converged || polish_converged`)
2509 /// would then surface that pixel as a "converged" fit in the spatial
2510 /// map. The guard at the top of `damped_fisher_stage` now keys
2511 /// convergence off `d_current.is_finite()`.
2512 ///
2513 /// Mirrors `lm.rs::test_all_fixed_params_nan_model` (issue #125.1),
2514 /// which exercises the equivalent guard in
2515 /// `levenberg_marquardt_with_mask`.
2516 #[test]
2517 fn test_joint_poisson_all_fixed_nan_transmission_does_not_converge() {
2518 struct NanModel {
2519 n_e: usize,
2520 }
2521 impl FitModel for NanModel {
2522 fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
2523 Ok(vec![f64::NAN; self.n_e])
2524 }
2525 }
2526
2527 let n_bins = 5;
2528 let o = vec![10.0; n_bins];
2529 let s = vec![5.0; n_bins];
2530 let model = NanModel { n_e: n_bins };
2531 let obj = JointPoissonObjective {
2532 model: &model,
2533 o: &o,
2534 s: &s,
2535 c: 1.0,
2536 active_mask: None,
2537 };
2538 let mut params = ParameterSet::new(vec![FitParameter::fixed("T", 0.5)]);
2539 let cfg = JointPoissonFitConfig::default();
2540
2541 let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
2542
2543 assert!(
2544 r.deviance.is_nan(),
2545 "expected NaN deviance from all-fixed NaN model; got {}",
2546 r.deviance
2547 );
2548 assert!(
2549 r.deviance_per_dof.is_nan(),
2550 "expected NaN deviance_per_dof; got {}",
2551 r.deviance_per_dof
2552 );
2553 assert!(
2554 !r.gn_converged,
2555 "all-fixed NaN deviance must not be reported as GN-converged",
2556 );
2557 assert_eq!(r.n_free, 0);
2558 assert_eq!(r.n_active, n_bins);
2559 // The damped-Fisher loop increments `iter` before the `n_free == 0`
2560 // branch hits `break`, so the all-fixed path always reports exactly
2561 // one iteration. Lock that in so future loop refactors don't
2562 // silently drift the iteration count.
2563 assert_eq!(
2564 r.gn_iterations, 1,
2565 "all-fixed branch should report exactly one iteration",
2566 );
2567 }
2568
2569 /// Companion to [`test_joint_poisson_all_fixed_nan_transmission_does_not_converge`]
2570 /// covering the polish-enabled path.
2571 ///
2572 /// `nelder_mead_minimize` asserts that `x0` is non-empty (see
2573 /// `nelder_mead.rs`), which used to panic when stage 2 was invoked with
2574 /// every parameter fixed. The polish entry-point now short-circuits on
2575 /// `free_indices().is_empty()`, so the call must return cleanly with
2576 /// `polish_converged == false` and the stage-1 NaN deviance preserved.
2577 /// Mirrors the pipeline configuration in `nereids-pipeline` where
2578 /// `with_counts_enable_polish(Some(true))` is set independently of
2579 /// whether the parameter set has any free entries.
2580 #[test]
2581 fn test_joint_poisson_all_fixed_nan_transmission_with_polish_does_not_panic() {
2582 struct NanModel {
2583 n_e: usize,
2584 }
2585 impl FitModel for NanModel {
2586 fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
2587 Ok(vec![f64::NAN; self.n_e])
2588 }
2589 }
2590
2591 let n_bins = 5;
2592 let o = vec![10.0; n_bins];
2593 let s = vec![5.0; n_bins];
2594 let model = NanModel { n_e: n_bins };
2595 let obj = JointPoissonObjective {
2596 model: &model,
2597 o: &o,
2598 s: &s,
2599 c: 1.0,
2600 active_mask: None,
2601 };
2602 let mut params = ParameterSet::new(vec![FitParameter::fixed("T", 0.5)]);
2603 let cfg = JointPoissonFitConfig {
2604 enable_polish: true,
2605 ..JointPoissonFitConfig::default()
2606 };
2607
2608 // Must not panic — the empty-x0 guard short-circuits stage 2.
2609 let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
2610
2611 assert!(
2612 r.deviance.is_nan(),
2613 "expected NaN deviance from all-fixed NaN model; got {}",
2614 r.deviance
2615 );
2616 assert!(
2617 !r.gn_converged,
2618 "all-fixed NaN deviance must not be reported as GN-converged",
2619 );
2620 assert!(
2621 !r.polish_converged,
2622 "polish stage must report not-converged when skipped on all-fixed params",
2623 );
2624 assert!(
2625 !r.polish_improved,
2626 "polish stage cannot have improved the deviance when it was skipped",
2627 );
2628 assert_eq!(
2629 r.polish_iterations, 0,
2630 "polish stage must report zero iterations when skipped",
2631 );
2632 assert_eq!(r.n_free, 0);
2633 assert_eq!(r.n_active, n_bins);
2634 assert_eq!(
2635 r.gn_iterations, 1,
2636 "all-fixed branch should report exactly one iteration",
2637 );
2638 }
2639
2640 /// Polish path with at least one **free** parameter must not report
2641 /// `polish_converged = true` when stage 1 ended on a non-finite
2642 /// deviance.
2643 ///
2644 /// Without the `best_d_stage1.is_finite()` short-circuit in the polish
2645 /// guard, Nelder-Mead would still run and return a finite `nm.fun`
2646 /// (its infeasible-point handler maps NaN evaluations to `+∞` and
2647 /// contracts away from them). The commit test `nm.fun < best_d_stage1`
2648 /// then reduces to `finite < NaN == false`, so the polish step is
2649 /// discarded — but `polish_converged` would inherit `nm.self_converged`
2650 /// regardless, leaking a spurious converged flag together with a NaN
2651 /// final deviance. Downstream pipeline code (`pipeline.rs`'s
2652 /// `gn_converged || polish_converged`) would then surface that fit as
2653 /// converged in the spatial map.
2654 ///
2655 /// Symmetric to the all-fixed NaN guard above: stage 2 refuses to run
2656 /// when there is no finite stage-1 deviance to refine.
2657 #[test]
2658 fn test_joint_poisson_polish_does_not_report_converged_when_stage1_nan() {
2659 struct NanModel {
2660 n_e: usize,
2661 }
2662 impl FitModel for NanModel {
2663 fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
2664 Ok(vec![f64::NAN; self.n_e])
2665 }
2666 }
2667
2668 let n_bins = 5;
2669 let o = vec![10.0; n_bins];
2670 let s = vec![5.0; n_bins];
2671 let model = NanModel { n_e: n_bins };
2672 let obj = JointPoissonObjective {
2673 model: &model,
2674 o: &o,
2675 s: &s,
2676 c: 1.0,
2677 active_mask: None,
2678 };
2679 // At least one FREE parameter so polish actually runs (unlike
2680 // `test_joint_poisson_all_fixed_nan_transmission_with_polish_does_not_panic`,
2681 // which exercises the empty-free-set short-circuit instead).
2682 let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
2683 let cfg = JointPoissonFitConfig {
2684 enable_polish: true,
2685 ..JointPoissonFitConfig::default()
2686 };
2687
2688 let r = joint_poisson_fit(&obj, &mut params, &cfg).unwrap();
2689
2690 assert!(
2691 r.deviance.is_nan(),
2692 "expected NaN deviance from NaN model; got {}",
2693 r.deviance
2694 );
2695 assert!(!r.gn_converged, "stage 1 cannot converge on NaN deviance",);
2696 assert!(
2697 !r.polish_converged,
2698 "stage 2 must not report converged when stage 1 ended non-finite",
2699 );
2700 assert!(
2701 !r.polish_improved,
2702 "polish cannot have improved a NaN starting deviance",
2703 );
2704 assert_eq!(
2705 r.polish_iterations, 0,
2706 "polish must not run when stage 1 is non-finite",
2707 );
2708 assert_eq!(r.n_free, 1);
2709 assert_eq!(r.n_active, n_bins);
2710 }
2711
2712 // ==================================================================
2713 // NaN-in-Jacobian during FD probes (Fisher info).
2714 //
2715 // The post-convergence Fisher / covariance path builds a Jacobian
2716 // via FD when the model has no analytical form. If the FD probe
2717 // straddles a region where the model returns NaN, the resulting
2718 // column is poisoned and the inverse Fisher inherits NaN entries.
2719 // The main LM loop's trial guard does not run here (it only checks
2720 // the trial step in the main optimisation loop).
2721 //
2722 // Per-cell skip: when the FD probe output is non-finite, leave the
2723 // entry at its zero default rather than dividing NaN by `actual_step`
2724 // (consistent with the "model-evaluation-failed" branch in
2725 // `compute_jacobian`).
2726 // ==================================================================
2727
2728 /// `fisher_information_fd` zeroes per-cell entries whose FD probe
2729 /// returned a non-finite model output, rather than baking NaN into
2730 /// the Fisher matrix (and from there into the inverse covariance).
2731 #[test]
2732 fn test_fisher_information_fd_skips_nan_probe() {
2733 // Model: T_i = θ_0 (constant). Returns NaN whenever
2734 // |θ_0 - 0.6| > 1e-3 — i.e. a NaN ring around the FD probe,
2735 // but a finite value at the base point.
2736 struct NanFdProbe;
2737 impl FitModel for NanFdProbe {
2738 fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
2739 let t = if (params[0] - 0.6).abs() > 1e-3 {
2740 f64::NAN
2741 } else {
2742 params[0]
2743 };
2744 Ok(vec![t; 3])
2745 }
2746 // No analytical_jacobian -> Fisher info must use FD fallback.
2747 }
2748 let model = NanFdProbe;
2749 let n = 3;
2750 let o = vec![10.0; n];
2751 let s = vec![5.0; n];
2752 let obj = JointPoissonObjective {
2753 model: &model,
2754 o: &o,
2755 s: &s,
2756 c: 1.0,
2757 active_mask: None,
2758 };
2759 let mut params = ParameterSet::new(vec![FitParameter::non_negative("T", 0.6)]);
2760 let info = obj
2761 .fisher_information_fd(&mut params, 1e-2)
2762 .expect("fisher_information_fd should not return Err on a finite base")
2763 .expect("fisher_information_fd should return Some(matrix)");
2764 // Every entry must be finite — column was skipped on NaN probe.
2765 for v in info.data.iter() {
2766 assert!(
2767 v.is_finite(),
2768 "fisher_information_fd produced non-finite entry: {v}"
2769 );
2770 }
2771 }
2772
2773 // ==================================================================
2774 // Per-element count validation propagates through `validate_inputs`.
2775 //
2776 // Round-2 added `validate_counts` only at the `joint_poisson_fit`
2777 // entry point. Direct callers of `deviance_from_transmission` /
2778 // `fisher_information_fd` / `profile_lambda_per_bin` (diagnostics
2779 // paths) bypassed that check, so a NaN in `o` would propagate
2780 // straight into the deviance sum via `NaN <= 0.0 == false` slipping
2781 // past `xlogy_ratio`'s zero-branch, and a negative count would be
2782 // silently swallowed as zero. Round-3 lifts the per-element check
2783 // into `validate_inputs`, which every public method already calls.
2784 // These tests run in release mode (no `debug_assert!`) and verify
2785 // the typed error reaches the caller.
2786 // ==================================================================
2787
2788 /// `deviance_from_transmission` must reject a NaN open-beam count
2789 /// with `InvalidConfig` rather than returning `Ok(NaN)` (or, worse,
2790 /// `Ok(finite)` if a future `xlogy_ratio` rewrite handled NaN by
2791 /// falling through to the zero branch). Regression for Round-3
2792 /// fix #1 — the inner `debug_assert!` is a no-op in release builds.
2793 #[test]
2794 fn test_deviance_from_transmission_rejects_non_finite_counts() {
2795 let n_bins = 4;
2796 let mut o = vec![10.0; n_bins];
2797 o[2] = f64::NAN;
2798 let s = vec![5.0; n_bins];
2799 let model = ConstModel { n_e: n_bins };
2800 let obj = JointPoissonObjective {
2801 model: &model,
2802 o: &o,
2803 s: &s,
2804 c: 1.0,
2805 active_mask: None,
2806 };
2807 let t = vec![0.5; n_bins];
2808 let err = obj.deviance_from_transmission(&t).unwrap_err();
2809 assert!(
2810 matches!(err, FittingError::InvalidConfig(ref msg) if msg.contains("open_beam_counts")),
2811 "expected InvalidConfig naming open_beam_counts; got {err:?}"
2812 );
2813
2814 // +inf likewise.
2815 let mut s_inf = vec![5.0; n_bins];
2816 s_inf[0] = f64::INFINITY;
2817 let obj_inf = JointPoissonObjective {
2818 model: &model,
2819 o: &vec![10.0; n_bins],
2820 s: &s_inf,
2821 c: 1.0,
2822 active_mask: None,
2823 };
2824 let err = obj_inf.deviance_from_transmission(&t).unwrap_err();
2825 assert!(
2826 matches!(err, FittingError::InvalidConfig(ref msg) if msg.contains("sample_counts")),
2827 "expected InvalidConfig naming sample_counts; got {err:?}"
2828 );
2829 }
2830
2831 /// `deviance_from_transmission` must reject a negative count with
2832 /// `InvalidConfig` rather than silently treating it as a zero-count
2833 /// bin (which `xlogy_ratio`'s `x <= 0.0` branch would do). Negatives
2834 /// indicate an upstream loader / TOF-subtraction bug; swallowing
2835 /// them as "no data" conceals the failure mode.
2836 #[test]
2837 fn test_deviance_from_transmission_rejects_negative_counts() {
2838 let n_bins = 3;
2839 let mut o = vec![10.0; n_bins];
2840 o[1] = -2.0;
2841 let s = vec![5.0; n_bins];
2842 let model = ConstModel { n_e: n_bins };
2843 let obj = JointPoissonObjective {
2844 model: &model,
2845 o: &o,
2846 s: &s,
2847 c: 1.0,
2848 active_mask: None,
2849 };
2850 let t = vec![0.5; n_bins];
2851 let err = obj.deviance_from_transmission(&t).unwrap_err();
2852 assert!(
2853 matches!(err, FittingError::InvalidConfig(ref msg) if msg.contains("open_beam_counts")),
2854 "expected InvalidConfig naming open_beam_counts; got {err:?}"
2855 );
2856 }
2857
2858 /// The reorientation also reaches `profile_lambda_per_bin` and
2859 /// `fisher_information_fd`: every public method that calls
2860 /// `validate_inputs` now picks up the per-element check.
2861 #[test]
2862 fn test_other_public_methods_reject_non_finite_counts() {
2863 let n_bins = 4;
2864 let mut s = vec![5.0; n_bins];
2865 s[3] = f64::NAN;
2866 let o = vec![10.0; n_bins];
2867 let model = ConstModel { n_e: n_bins };
2868 let obj = JointPoissonObjective {
2869 model: &model,
2870 o: &o,
2871 s: &s,
2872 c: 1.0,
2873 active_mask: None,
2874 };
2875 let t = vec![0.5; n_bins];
2876
2877 let err = obj.profile_lambda_per_bin(&t).unwrap_err();
2878 assert!(
2879 matches!(err, FittingError::InvalidConfig(_)),
2880 "profile_lambda_per_bin: expected InvalidConfig; got {err:?}"
2881 );
2882
2883 let params = vec![0.5];
2884 let free_idx = vec![0];
2885 let err = obj
2886 .deviance_gradient_analytical(¶ms, &free_idx)
2887 .unwrap_err();
2888 assert!(
2889 matches!(err, FittingError::InvalidConfig(_)),
2890 "deviance_gradient_analytical: expected InvalidConfig; got {err:?}"
2891 );
2892
2893 let err = obj.fisher_information(¶ms, &free_idx).unwrap_err();
2894 assert!(
2895 matches!(err, FittingError::InvalidConfig(_)),
2896 "fisher_information: expected InvalidConfig; got {err:?}"
2897 );
2898
2899 let mut ps = ParameterSet::new(vec![FitParameter::non_negative("T", 0.5)]);
2900 let err = obj.fisher_information_fd(&mut ps, 1e-2).unwrap_err();
2901 assert!(
2902 matches!(err, FittingError::InvalidConfig(_)),
2903 "fisher_information_fd: expected InvalidConfig; got {err:?}"
2904 );
2905 }
2906
2907 /// `validate_inputs` now reports caller-supplied transmission length
2908 /// mismatches with `field = "transmission"` and `expected = o.len()`.
2909 /// Pre-fix this used `field = "open_beam_counts"` with reversed
2910 /// expected/actual, which read as "the open-beam array is wrong"
2911 /// when the actual fault was the caller's `t` slice. Regression
2912 /// for Round-3 fix #2.
2913 #[test]
2914 fn test_validate_inputs_reports_transmission_length_mismatch_correctly() {
2915 let n_bins = 5;
2916 let o = vec![10.0; n_bins];
2917 let s = vec![5.0; n_bins];
2918 let model = ConstModel { n_e: n_bins };
2919 let obj = JointPoissonObjective {
2920 model: &model,
2921 o: &o,
2922 s: &s,
2923 c: 1.0,
2924 active_mask: None,
2925 };
2926 // Caller passes `t` shorter than `o`/`s`.
2927 let t_short = vec![0.5; n_bins - 2];
2928 let err = obj.deviance_from_transmission(&t_short).unwrap_err();
2929 match err {
2930 FittingError::LengthMismatch {
2931 expected,
2932 actual,
2933 field,
2934 } => {
2935 assert_eq!(field, "transmission", "field must name `transmission`");
2936 assert_eq!(expected, n_bins, "expected must be o.len()");
2937 assert_eq!(actual, n_bins - 2, "actual must be t.len()");
2938 }
2939 other => panic!("expected LengthMismatch on transmission; got {other:?}"),
2940 }
2941 }
2942}