Skip to main content

nereids_fitting/
lm.rs

1//! Levenberg-Marquardt least-squares optimizer.
2//!
3//! Minimizes χ² = Σᵢ [(y_obs - y_model)² / σᵢ²] by iteratively solving:
4//!
5//!   (JᵀWJ + λ·diag(JᵀWJ)) · δ = JᵀW·r
6//!
7//! where J is the Jacobian, W = diag(1/σ²), r = y_obs - y_model,
8//! and λ is the damping parameter.
9//!
10//! ## SAMMY Reference
11//! - `fit/` module, manual Sec. IV (Bayes equations / generalized least-squares)
12
13use nereids_core::constants::{LM_DIAGONAL_FLOOR, PIVOT_FLOOR};
14
15use crate::error::FittingError;
16use crate::parameters::ParameterSet;
17
18/// Row-major flat matrix for cache-friendly storage.
19///
20/// Replaces `Vec<Vec<f64>>` to collapse ~N separate heap allocations into 1
21/// and improve cache locality for JtWJ assembly.  Access: `data[i * ncols + j]`.
22#[derive(Debug, Clone)]
23pub struct FlatMatrix {
24    /// Flat row-major storage: `data[i * ncols + j]` = element at row i, col j.
25    pub data: Vec<f64>,
26    /// Number of rows.
27    pub nrows: usize,
28    /// Number of columns.
29    pub ncols: usize,
30}
31
32impl FlatMatrix {
33    /// Create a new zero-filled matrix with the given dimensions.
34    pub fn zeros(nrows: usize, ncols: usize) -> Self {
35        let len = nrows
36            .checked_mul(ncols)
37            .expect("FlatMatrix dimensions overflow usize");
38        Self {
39            data: vec![0.0; len],
40            nrows,
41            ncols,
42        }
43    }
44
45    /// Access element at (row, col) immutably.
46    #[inline(always)]
47    pub fn get(&self, row: usize, col: usize) -> f64 {
48        debug_assert!(row < self.nrows && col < self.ncols);
49        self.data[row * self.ncols + col]
50    }
51
52    /// Access element at (row, col) mutably.
53    #[inline(always)]
54    pub fn get_mut(&mut self, row: usize, col: usize) -> &mut f64 {
55        debug_assert!(row < self.nrows && col < self.ncols);
56        &mut self.data[row * self.ncols + col]
57    }
58}
59
60/// #125.4: Maximum damping parameter before the optimizer gives up.
61///
62/// When λ exceeds this threshold, the optimizer is stuck in a region where no
63/// step improves chi-squared.  Breaking out avoids wasting iterations.
64const LAMBDA_BREAKOUT: f64 = 1e16;
65
66/// Configuration for the LM optimizer.
67#[derive(Debug, Clone)]
68pub struct LmConfig {
69    /// Maximum number of iterations.
70    pub max_iter: usize,
71    /// Initial damping parameter λ.
72    pub lambda_init: f64,
73    /// Factor to increase λ on rejected step.
74    pub lambda_up: f64,
75    /// Factor to decrease λ on accepted step.
76    pub lambda_down: f64,
77    /// Convergence tolerance on relative χ² change.
78    pub tol_chi2: f64,
79    /// Convergence tolerance on relative parameter change.
80    pub tol_param: f64,
81    /// Step size for finite-difference Jacobian.
82    pub fd_step: f64,
83    /// Whether to compute the covariance matrix (and uncertainties) after
84    /// convergence.  This requires an extra Jacobian evaluation + matrix
85    /// inversion at the final parameters.  Set to `false` for per-pixel
86    /// spatial mapping where only densities are needed.
87    ///
88    /// Default: `true`.
89    pub compute_covariance: bool,
90}
91
92impl Default for LmConfig {
93    fn default() -> Self {
94        Self {
95            max_iter: 200,
96            lambda_init: 1e-3,
97            lambda_up: 10.0,
98            lambda_down: 0.1,
99            tol_chi2: 1e-8,
100            tol_param: 1e-8,
101            fd_step: 1e-6,
102            compute_covariance: true,
103        }
104    }
105}
106
107/// Result of a Levenberg-Marquardt fit.
108#[derive(Debug, Clone)]
109pub struct LmResult {
110    /// Final chi-squared value.
111    pub chi_squared: f64,
112    /// Reduced chi-squared (χ²/ν where ν = n_data - n_params).
113    pub reduced_chi_squared: f64,
114    /// Number of iterations taken.
115    pub iterations: usize,
116    /// Whether the fit converged.
117    pub converged: bool,
118    /// Final parameter values (all parameters, including fixed).
119    pub params: Vec<f64>,
120    /// Covariance matrix of free parameters (n_free × n_free), if available.
121    pub covariance: Option<FlatMatrix>,
122    /// Standard errors of free parameters (diagonal of covariance).
123    pub uncertainties: Option<Vec<f64>>,
124}
125
126/// A model function that can be fitted.
127///
128/// Given parameter values (all params including fixed), computes
129/// the model prediction at each data point.
130pub trait FitModel {
131    /// Evaluate the model for the given parameters.
132    ///
133    /// On success, returns a vector of model predictions with the same
134    /// length as the data being fitted. On failure, returns a
135    /// [`FittingError`] indicating that the model could not be evaluated
136    /// (e.g. a broadening or physics error).
137    ///
138    /// **Optimizer semantics:** during Levenberg-Marquardt trial steps,
139    /// an `Err` is treated as a failed step (increase λ / backtrack).
140    /// At the initial point or post-convergence, `Err` is propagated.
141    fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError>;
142
143    /// Optionally provide an analytical Jacobian.
144    ///
145    /// `free_param_indices`: indices (into `params`) of the free parameters,
146    /// in the same order as the Jacobian columns.
147    ///
148    /// `y_current`: current model output, i.e. `self.evaluate(params)`.
149    /// Provided so implementations can compute J analytically from T without
150    /// an extra `evaluate` call.
151    ///
152    /// Returns `Some(J)` where `J.get(i, j) = ∂model[i]/∂params[free_param_indices[j]]`.
153    /// The matrix has `y_current.len()` rows and `free_param_indices.len()` columns.
154    /// Return `None` to fall back to finite-difference Jacobian (the default).
155    fn analytical_jacobian(
156        &self,
157        _params: &[f64],
158        _free_param_indices: &[usize],
159        _y_current: &[f64],
160    ) -> Option<FlatMatrix> {
161        None
162    }
163}
164
165/// Blanket implementation: shared references to any `FitModel` also implement
166/// `FitModel`, forwarding all calls to the underlying implementation.
167///
168/// This enables `NormalizedTransmissionModel<&dyn FitModel>` to work when the
169/// inner model is borrowed (e.g. in `fit_spectrum` where the inner model is a
170/// local variable wrapped conditionally).
171impl<M: FitModel + ?Sized> FitModel for &M {
172    fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
173        (**self).evaluate(params)
174    }
175
176    fn analytical_jacobian(
177        &self,
178        params: &[f64],
179        free_param_indices: &[usize],
180        y_current: &[f64],
181    ) -> Option<FlatMatrix> {
182        (**self).analytical_jacobian(params, free_param_indices, y_current)
183    }
184}
185
186/// Compute weighted chi-squared: Σ [(y_obs - y_model)² / σ²].
187///
188/// When `active_mask` is `Some(m)`, bins where `m[i]` is `false` are
189/// **skipped entirely** rather than zero-weighted.  Explicit row-skip
190/// matters when masked bins may contain non-finite residuals (e.g. NaN
191/// outside the user's fit-energy range): `0.0 * NaN = NaN`, so a
192/// zero-weight strategy would propagate NaN through the χ² accumulator
193/// even though the bin contributes no information.  See SAMMY EMIN/EMAX
194/// semantics (#514) and lm.rs Fix 4.
195fn chi_squared(residuals: &[f64], weights: &[f64], active_mask: Option<&[bool]>) -> f64 {
196    let mut sum = 0.0;
197    for (i, (&r, &w)) in residuals.iter().zip(weights.iter()).enumerate() {
198        if active_mask.is_some_and(|m| !m[i]) {
199            continue;
200        }
201        sum += r * r * w;
202    }
203    sum
204}
205
206/// Infinity norm of the gradient, scaled by local curvature and residual size.
207///
208/// This is a dimensionless first-order optimality measure for least squares:
209/// small values indicate that the current point is stationary even when χ² is
210/// nonzero because the data are noisy or the model is imperfect.
211fn scaled_gradient_inf_norm(jtw_j: &FlatMatrix, jtw_r: &[f64], chi2: f64) -> f64 {
212    let residual_scale = chi2.sqrt().max(1.0);
213    let mut max_scaled: f64 = 0.0;
214    for (j, &grad_j) in jtw_r.iter().enumerate() {
215        let curvature = jtw_j.get(j, j).abs().sqrt();
216        let scale = curvature * residual_scale + PIVOT_FLOOR;
217        max_scaled = max_scaled.max(grad_j.abs() / scale);
218    }
219    max_scaled
220}
221
222/// Whether the local model has meaningful curvature in at least one free direction.
223///
224/// Purely flat models (J = 0 everywhere) should still report failure rather than
225/// "converged", even though their gradient is numerically zero.
226fn has_informative_curvature(jtw_j: &FlatMatrix) -> bool {
227    (0..jtw_j.nrows).any(|j| jtw_j.get(j, j) > LM_DIAGONAL_FLOOR)
228}
229
230/// Compute the Jacobian, preferring an analytical formula over finite differences.
231///
232/// `y_current` must equal `model.evaluate(&params.all_values())` at the
233/// current parameter values — it is passed in to avoid a redundant evaluate
234/// call (the LM loop already has this vector from the previous accepted step).
235///
236/// `all_vals_buf` is a scratch buffer reused across the per-parameter FD loop
237/// to avoid allocating a fresh `Vec<f64>` on every `model.evaluate()` call.
238///
239/// `free_idx_buf` is a scratch buffer for `params.free_indices_into()`, reused
240/// across iterations to avoid per-Jacobian allocation.
241///
242/// J.get(i, j) = ∂model[i] / ∂free_param[j]
243fn compute_jacobian(
244    model: &dyn FitModel,
245    params: &mut ParameterSet,
246    y_current: &[f64],
247    fd_step: f64,
248    all_vals_buf: &mut Vec<f64>,
249    free_idx_buf: &mut Vec<usize>,
250) -> Result<FlatMatrix, FittingError> {
251    params.free_indices_into(free_idx_buf);
252    let n_free = free_idx_buf.len();
253    let n_data = y_current.len();
254
255    // Try analytical Jacobian first (no extra evaluate calls).
256    params.all_values_into(all_vals_buf);
257    if let Some(j) = model.analytical_jacobian(all_vals_buf, free_idx_buf, y_current) {
258        debug_assert!(
259            j.nrows == n_data && j.ncols == n_free && j.data.len() == n_data * n_free,
260            "analytical_jacobian shape mismatch: got ({}x{}, len={}), expected ({}x{}, len={})",
261            j.nrows,
262            j.ncols,
263            j.data.len(),
264            n_data,
265            n_free,
266            n_data * n_free,
267        );
268        return Ok(j);
269    }
270
271    // Fallback: forward finite differences, reusing y_current as the base.
272    let mut jacobian = FlatMatrix::zeros(n_data, n_free);
273
274    for (j, &idx) in free_idx_buf.iter().enumerate() {
275        let original = params.params[idx].value;
276        let step = fd_step * (1.0 + original.abs());
277
278        params.params[idx].value = original + step;
279        params.params[idx].clamp();
280        let mut actual_step = params.params[idx].value - original;
281
282        // #112: If the forward step is blocked by an upper bound, try the
283        // backward step so the Jacobian column is not frozen at zero.
284        if actual_step.abs() < PIVOT_FLOOR {
285            params.params[idx].value = original - step;
286            params.params[idx].clamp();
287            actual_step = params.params[idx].value - original;
288            if actual_step.abs() < PIVOT_FLOOR {
289                // Truly stuck at a point constraint — skip this parameter.
290                params.params[idx].value = original;
291                continue;
292            }
293        }
294
295        params.all_values_into(all_vals_buf);
296        let perturbed = match model.evaluate(all_vals_buf) {
297            Ok(v) => v,
298            Err(_) => {
299                // Restore original before skipping — leaving the column as
300                // zero makes this parameter unresponsive for one LM step,
301                // which is safe (matches Poisson compute_gradient pattern).
302                params.params[idx].value = original;
303                continue;
304            }
305        };
306        params.params[idx].value = original;
307
308        // The main LM loop checks the trial step via
309        // `trial_has_active_nonfinite`, but `compute_jacobian` is also
310        // called from the post-convergence covariance path, where no
311        // such guard runs.  A NaN row in the perturbed output divided
312        // by `actual_step` would yield NaN in the Jacobian, which then
313        // poisons JᵀWJ and the inverse covariance.
314        //
315        // Per-cell skip (zero the entry) rather than whole-column skip:
316        // a NaN at a masked / inactive row is benign (the JᵀWJ assembly
317        // skips masked rows entirely via `active_mask.is_some_and(|m| !m[i])`),
318        // so a finite Jacobian row produced from a bad-but-masked probe
319        // must still be allowed to land in the column — see
320        // `test_lm_active_mask_tolerates_model_nan_outside_range`.  The
321        // active-row NaNs that would otherwise propagate into JᵀWJ are
322        // zeroed here, which is the same numerical outcome as if the
323        // model had reported `Err` at that probe (skipped column, but
324        // per-row).
325        for i in 0..n_data {
326            let p = perturbed[i];
327            let y = y_current[i];
328            if p.is_finite() && y.is_finite() {
329                *jacobian.get_mut(i, j) = (p - y) / actual_step;
330            }
331            // else: leave at the zero-default; for active rows this is
332            // safe (matches the whole-column-skip outcome on Err), for
333            // masked rows it never gets read.
334        }
335    }
336
337    Ok(jacobian)
338}
339
340/// Solve (A + λ·diag(A)) · x = b using Gaussian elimination.
341///
342/// A is a flat n×n symmetric positive definite matrix (approximately).
343/// Returns the solution vector x.
344pub(crate) fn solve_damped_system(a: &FlatMatrix, b: &[f64], lambda: f64) -> Option<Vec<f64>> {
345    let n = b.len();
346    if n == 0 {
347        return Some(vec![]);
348    }
349
350    // Build the augmented matrix [A + λ·diag(A) | b] as flat (n × (n+1)).
351    let ncols = n + 1;
352    let mut aug = FlatMatrix::zeros(n, ncols);
353    for (i, &bi) in b.iter().enumerate() {
354        for j in 0..n {
355            *aug.get_mut(i, j) = a.get(i, j);
356        }
357        *aug.get_mut(i, i) += lambda * a.get(i, i).max(LM_DIAGONAL_FLOOR); // Ensure non-zero diagonal
358        *aug.get_mut(i, n) = bi;
359    }
360
361    // Gaussian elimination with partial pivoting
362    for col in 0..n {
363        // Find pivot
364        let mut max_val = aug.get(col, col).abs();
365        let mut max_row = col;
366        for row in (col + 1)..n {
367            if aug.get(row, col).abs() > max_val {
368                max_val = aug.get(row, col).abs();
369                max_row = row;
370            }
371        }
372
373        if max_val < PIVOT_FLOOR {
374            return None; // Singular
375        }
376
377        // Swap rows col and max_row in the flat buffer.
378        if col != max_row {
379            let (row_a, row_b) = (col * ncols, max_row * ncols);
380            let (first, second) = aug.data.split_at_mut(row_b);
381            first[row_a..row_a + ncols].swap_with_slice(&mut second[..ncols]);
382        }
383
384        let pivot = aug.get(col, col);
385        for row in (col + 1)..n {
386            let factor = aug.get(row, col) / pivot;
387            for j in col..=n {
388                let val = aug.get(col, j);
389                *aug.get_mut(row, j) -= factor * val;
390            }
391        }
392    }
393
394    // Back substitution
395    let mut x = vec![0.0; n];
396    for i in (0..n).rev() {
397        let mut sum = aug.get(i, n);
398        for (j, &xj) in x.iter().enumerate().skip(i + 1) {
399            sum -= aug.get(i, j) * xj;
400        }
401        x[i] = sum / aug.get(i, i);
402    }
403
404    Some(x)
405}
406
407/// Invert a symmetric positive definite matrix (for covariance).
408///
409/// Input: flat n×n matrix. Output: flat n×n inverse, or None if singular.
410pub(crate) fn invert_matrix(a: &FlatMatrix) -> Option<FlatMatrix> {
411    let n = a.nrows;
412    if n == 0 {
413        return Some(FlatMatrix::zeros(0, 0));
414    }
415
416    // Build [A | I] as flat (n × 2n).
417    let ncols = 2 * n;
418    let mut aug = FlatMatrix::zeros(n, ncols);
419    for i in 0..n {
420        for j in 0..n {
421            *aug.get_mut(i, j) = a.get(i, j);
422        }
423        *aug.get_mut(i, n + i) = 1.0;
424    }
425
426    // Forward elimination with partial pivoting
427    for col in 0..n {
428        let mut max_val = aug.get(col, col).abs();
429        let mut max_row = col;
430        for row in (col + 1)..n {
431            if aug.get(row, col).abs() > max_val {
432                max_val = aug.get(row, col).abs();
433                max_row = row;
434            }
435        }
436
437        if max_val < PIVOT_FLOOR {
438            return None;
439        }
440
441        // Swap rows col and max_row.
442        if col != max_row {
443            let (row_a, row_b) = (col * ncols, max_row * ncols);
444            let (first, second) = aug.data.split_at_mut(row_b);
445            first[row_a..row_a + ncols].swap_with_slice(&mut second[..ncols]);
446        }
447
448        let pivot = aug.get(col, col);
449        for j in 0..ncols {
450            *aug.get_mut(col, j) /= pivot;
451        }
452
453        for row in 0..n {
454            if row != col {
455                let factor = aug.get(row, col);
456                for j in 0..ncols {
457                    let val = aug.get(col, j);
458                    *aug.get_mut(row, j) -= factor * val;
459                }
460            }
461        }
462    }
463
464    // Extract the right half [I|A⁻¹] → A⁻¹
465    let mut inv = FlatMatrix::zeros(n, n);
466    for i in 0..n {
467        for j in 0..n {
468            *inv.get_mut(i, j) = aug.get(i, n + j);
469        }
470    }
471
472    Some(inv)
473}
474
475/// Run the Levenberg-Marquardt optimizer.
476///
477/// # Arguments
478/// * `model` — Forward model implementing `FitModel`.
479/// * `y_obs` — Observed data values.
480/// * `sigma` — Uncertainties on observed data (standard deviations).
481/// * `params` — Initial parameter set (modified in place on convergence).
482/// * `config` — LM configuration.
483///
484/// # Returns
485/// Fit result including final parameters, chi-squared, and uncertainties.
486pub fn levenberg_marquardt(
487    model: &dyn FitModel,
488    y_obs: &[f64],
489    sigma: &[f64],
490    params: &mut ParameterSet,
491    config: &LmConfig,
492) -> Result<LmResult, FittingError> {
493    levenberg_marquardt_with_mask(model, y_obs, sigma, params, config, None)
494}
495
496/// Run the Levenberg-Marquardt optimizer with an optional per-bin
497/// active mask (SAMMY EMIN/EMAX-equivalent fit-energy-range restriction).
498///
499/// When `active_mask` is `Some(m)`:
500/// - bins where `m[i]` is `false` are excluded from χ² and the normal
501///   equations (weight zeroed before assembly), so they contribute
502///   nothing to the gradient or Hessian;
503/// - the model is still evaluated on the full grid so resolution
504///   broadening at the boundaries of the active region is correct;
505/// - `dof = n_active − n_free` where `n_active` is the count of
506///   active bins.
507///
508/// When `active_mask` is `None`, behaviour is identical to
509/// [`levenberg_marquardt`] (all bins active).
510pub fn levenberg_marquardt_with_mask(
511    model: &dyn FitModel,
512    y_obs: &[f64],
513    sigma: &[f64],
514    params: &mut ParameterSet,
515    config: &LmConfig,
516    active_mask: Option<&[bool]>,
517) -> Result<LmResult, FittingError> {
518    let n_data = y_obs.len();
519    if n_data == 0 {
520        return Err(FittingError::EmptyData);
521    }
522    if sigma.len() != n_data {
523        return Err(FittingError::LengthMismatch {
524            expected: n_data,
525            actual: sigma.len(),
526            field: "sigma",
527        });
528    }
529    if let Some(m) = active_mask
530        && m.len() != n_data
531    {
532        return Err(FittingError::LengthMismatch {
533            expected: n_data,
534            actual: m.len(),
535            field: "active_mask",
536        });
537    }
538    let n_active = crate::active_mask::active_count(active_mask, n_data);
539
540    // SAMMY EMIN/EMAX-equivalent fit-energy-range (#514): a mask with
541    // zero active bins means the user's `[E_min, E_max]` does not
542    // overlap the energy grid.  No data contributes to the cost
543    // function — return non-converged with NaN χ² rather than
544    // falling through to the all-fixed fast-return path (which would
545    // report `converged: true, chi_squared: 0` from a zero-row sum)
546    // or to the main optimisation loop (where `n_active < n_free`
547    // would catch it but only after wasted setup work).
548    if n_active == 0 {
549        return Ok(LmResult {
550            chi_squared: f64::NAN,
551            reduced_chi_squared: f64::NAN,
552            iterations: 0,
553            converged: false,
554            params: params.all_values(),
555            covariance: None,
556            uncertainties: None,
557        });
558    }
559
560    let n_free = params.n_free();
561
562    // Early return when all parameters are fixed: evaluate once and report the
563    // model's chi-squared.  There is nothing to optimize, so iterating would
564    // waste cycles.
565    if n_free == 0 {
566        let weights: Vec<f64> = sigma
567            .iter()
568            .enumerate()
569            .map(|(i, &s)| {
570                // Active-bin masking (SAMMY EMIN/EMAX): bins outside the
571                // user range contribute zero to χ².
572                if active_mask.is_some_and(|m| !m[i]) {
573                    return 0.0;
574                }
575                if !s.is_finite() || s <= 0.0 {
576                    1.0 / 1e30
577                } else {
578                    1.0 / (s * s)
579                }
580            })
581            .collect();
582        let y_model = model.evaluate(&params.all_values())?;
583
584        // #P1: If the model produces NaN/Inf with all-fixed parameters,
585        // return converged=false rather than silently propagating NaN chi².
586        // Covariance/uncertainties are None because the fit did not converge —
587        // an unconverged result has no meaningful covariance to report.
588        //
589        // SAMMY EMIN/EMAX-equivalent fit-energy-range (#514): masked rows are
590        // skipped throughout the cost / normal-equation accumulators, so a
591        // non-finite model output at a masked bin should not abort the fit
592        // — only check finiteness on active rows.
593        let model_has_active_nonfinite = y_model
594            .iter()
595            .enumerate()
596            .any(|(i, v)| active_mask.is_none_or(|m| m[i]) && !v.is_finite());
597        if model_has_active_nonfinite {
598            return Ok(LmResult {
599                chi_squared: f64::NAN,
600                reduced_chi_squared: f64::NAN,
601                iterations: 0,
602                converged: false,
603                params: params.all_values(),
604                covariance: None,
605                uncertainties: None,
606            });
607        }
608
609        let residuals: Vec<f64> = y_obs
610            .iter()
611            .zip(y_model.iter())
612            .map(|(&obs, &mdl)| obs - mdl)
613            .collect();
614        let chi2 = chi_squared(&residuals, &weights, active_mask);
615        // #125.5: Compute dof via `n_active - n_free` (with `n_active = n_data`
616        // when no mask is set) to mirror the main path and keep a single
617        // visible formula.  When a fit-energy-range mask is in effect,
618        // `n_active` reflects the count of bins that actually contributed
619        // to χ² (SAMMY EMIN/EMAX semantics).
620        let dof = n_active.saturating_sub(n_free);
621        let reduced = if dof > 0 { chi2 / dof as f64 } else { f64::NAN };
622        return Ok(LmResult {
623            chi_squared: chi2,
624            reduced_chi_squared: reduced,
625            iterations: 0,
626            converged: true,
627            params: params.all_values(),
628            covariance: Some(FlatMatrix::zeros(0, 0)),
629            uncertainties: Some(vec![]),
630        });
631    }
632
633    // #108.3: Underdetermined systems — when n_active < n_free, the problem is
634    // underdetermined and the Jacobian cannot be full rank.  Return early
635    // with converged=false so callers can detect the problem.
636    // n_active == n_free is exactly determined (dof=0) and still solvable;
637    // we allow it and report reduced_chi_squared = NaN (0/0).
638    //
639    // When a fit-energy-range mask is in effect, the underdetermined
640    // check uses the active-bin count (SAMMY EMIN/EMAX semantics) — masked
641    // bins do not contribute information to the fit.
642    if n_active < n_free {
643        return Ok(LmResult {
644            chi_squared: f64::NAN,
645            reduced_chi_squared: f64::NAN,
646            iterations: 0,
647            converged: false,
648            params: params.all_values(),
649            covariance: None,
650            uncertainties: None,
651        });
652    }
653    let dof = n_active - n_free;
654
655    // #104: Validate sigma — division by zero or non-finite sigma would produce
656    // NaN/Inf weights and silently corrupt the entire fit.  Clamp to a small
657    // floor instead of rejecting outright, so callers with a few zero-sigma
658    // bins still get a usable fit.
659    //
660    // Active-bin masking (SAMMY EMIN/EMAX): bins outside the user range
661    // get weight 0 here, which propagates through χ² and the JᵀWJ /
662    // JᵀWr normal-equation assembly to contribute exactly zero —
663    // covering both residual and gradient masking with no further
664    // changes downstream.
665    let weights: Vec<f64> = sigma
666        .iter()
667        .enumerate()
668        .map(|(i, &s)| {
669            if active_mask.is_some_and(|m| !m[i]) {
670                return 0.0;
671            }
672            if !s.is_finite() || s <= 0.0 {
673                // Treat as negligible weight (huge sigma) rather than panicking.
674                1.0 / 1e30
675            } else {
676                1.0 / (s * s)
677            }
678        })
679        .collect();
680
681    // Scratch buffers reused across the optimization loop for
682    // params.all_values_into() calls in compute_jacobian (1 + N_free calls
683    // per Jacobian computation) and the trial-step evaluation,
684    // params.free_values_into() calls for snapshotting free parameters
685    // before trial steps, and params.free_indices_into() calls inside
686    // compute_jacobian.
687    let mut all_vals_buf = Vec::with_capacity(params.params.len());
688    let mut free_vals_buf = Vec::with_capacity(n_free);
689    let mut free_idx_buf = Vec::with_capacity(n_free);
690
691    // Initial model output, residuals, and chi².
692    // y_current is kept up-to-date after accepted steps so that the next
693    // Jacobian call can reuse it without an extra evaluate() call.
694    params.all_values_into(&mut all_vals_buf);
695    let mut y_current = model.evaluate(&all_vals_buf)?;
696    let mut residuals: Vec<f64> = y_obs
697        .iter()
698        .zip(y_current.iter())
699        .map(|(&obs, &mdl)| obs - mdl)
700        .collect();
701    let mut chi2 = chi_squared(&residuals, &weights, active_mask);
702
703    let mut lambda = config.lambda_init;
704    let mut converged = false;
705    let mut iter = 0;
706
707    for _ in 0..config.max_iter {
708        iter += 1;
709
710        // Compute Jacobian — uses y_current to avoid a redundant evaluate().
711        // Analytical Jacobian (if provided by the model) costs 0 extra evaluates;
712        // finite-difference fallback costs N_free extra evaluates.
713        let jacobian = compute_jacobian(
714            model,
715            params,
716            &y_current,
717            config.fd_step,
718            &mut all_vals_buf,
719            &mut free_idx_buf,
720        )?;
721
722        // Build normal equations: JᵀWJ and JᵀWr.
723        //
724        // Active-bin masking (SAMMY EMIN/EMAX, #514): explicitly skip
725        // masked rows rather than relying on weights[i] == 0 — the
726        // model or the residual at a masked margin bin may be NaN
727        // (e.g. when y_obs is NaN outside the user's fit-energy range),
728        // and `0.0 * NaN = NaN` would poison both accumulators.
729        let mut jtw_j = FlatMatrix::zeros(n_free, n_free);
730        let mut jtw_r = vec![0.0; n_free];
731
732        for (i, (&wi, &ri)) in weights.iter().zip(residuals.iter()).enumerate() {
733            if active_mask.is_some_and(|m| !m[i]) {
734                continue;
735            }
736            for (j, jtw_r_j) in jtw_r.iter_mut().enumerate() {
737                let jij = jacobian.get(i, j);
738                *jtw_r_j += jij * wi * ri;
739                for k in 0..n_free {
740                    *jtw_j.get_mut(j, k) += jij * wi * jacobian.get(i, k);
741                }
742            }
743        }
744        let scaled_grad_inf = scaled_gradient_inf_norm(&jtw_j, &jtw_r, chi2);
745        let informative_curvature = has_informative_curvature(&jtw_j);
746
747        // Solve (JᵀWJ + λ·diag(JᵀWJ)) · δ = JᵀWr
748        let delta = match solve_damped_system(&jtw_j, &jtw_r, lambda) {
749            Some(d) => d,
750            None => break, // Singular system
751        };
752
753        // Trial step — snapshot free values into a reusable buffer to avoid
754        // per-iteration allocation.
755        params.free_values_into(&mut free_vals_buf);
756        let trial_free: Vec<f64> = free_vals_buf
757            .iter()
758            .zip(delta.iter())
759            .map(|(&v, &d)| v + d)
760            .collect();
761        let param_change: f64 = delta
762            .iter()
763            .zip(free_vals_buf.iter())
764            .map(|(&d, &v)| (d / (v.abs() + PIVOT_FLOOR)).powi(2))
765            .sum::<f64>()
766            .sqrt();
767        params.set_free_values(&trial_free);
768
769        params.all_values_into(&mut all_vals_buf);
770        let y_trial = match model.evaluate(&all_vals_buf) {
771            Ok(y) => y,
772            Err(_) => {
773                // Treat evaluation error as a bad step (same as NaN).
774                params.set_free_values(&free_vals_buf);
775                lambda *= config.lambda_up;
776                if lambda > LAMBDA_BREAKOUT {
777                    converged = false;
778                    break;
779                }
780                continue;
781            }
782        };
783
784        // #113: If the model produced NaN/Inf at any *active* bin, treat as
785        // a bad step (same as chi2 increase) — increase lambda and try again.
786        // SAMMY EMIN/EMAX-equivalent fit-energy-range (#514): masked rows are
787        // skipped from χ² / JᵀWJ / JᵀWr / covariance, so a non-finite model
788        // output at a masked margin bin must not penalise the trial step.
789        let trial_has_active_nonfinite = y_trial
790            .iter()
791            .enumerate()
792            .any(|(i, v)| active_mask.is_none_or(|m| m[i]) && !v.is_finite());
793        if trial_has_active_nonfinite {
794            params.set_free_values(&free_vals_buf);
795            lambda *= config.lambda_up;
796            if lambda > LAMBDA_BREAKOUT {
797                converged = false;
798                break;
799            }
800            continue;
801        }
802
803        let trial_residuals: Vec<f64> = y_obs
804            .iter()
805            .zip(y_trial.iter())
806            .map(|(&obs, &mdl)| obs - mdl)
807            .collect();
808        let trial_chi2 = chi_squared(&trial_residuals, &weights, active_mask);
809        let chi2_delta = (trial_chi2 - chi2).abs();
810        let chi2_scale = chi2.abs().max(trial_chi2.abs()).max(1.0);
811        let chi2_stagnated = chi2_delta <= config.tol_chi2 * chi2_scale;
812
813        if trial_chi2 < chi2 {
814            // Accept step — cache y_trial so the next iteration can skip
815            // the base evaluate() inside compute_jacobian.
816            let rel_change = (chi2 - trial_chi2) / (chi2 + PIVOT_FLOOR);
817            chi2 = trial_chi2;
818            residuals = trial_residuals;
819            y_current = y_trial;
820            lambda *= config.lambda_down;
821
822            // Check convergence: relative chi2 change is tiny or parameters
823            // stopped moving.  The old third condition
824            // `chi2 < tol_chi2 * n_data` was scale-dependent and could cause
825            // premature convergence on data with small residuals.  (#108.2)
826            if rel_change < config.tol_chi2 || param_change < config.tol_param {
827                converged = true;
828                break;
829            }
830        } else {
831            // Numerical stagnation: the strict LM acceptance test keeps
832            // `trial_chi2 == chi2` in the reject path, but when both the
833            // objective change and parameter step are tiny, the optimizer may
834            // already be at a nonzero-χ² stationary point (noisy data,
835            // correlated parameters, imperfect model). Report convergence if
836            // the gradient is also tiny and the local model has real
837            // curvature, rather than inflating lambda until breakout.
838            let grad_tol = config.tol_chi2.sqrt().max(config.tol_param.sqrt());
839            if chi2_stagnated
840                && param_change < config.tol_param
841                && scaled_grad_inf < grad_tol
842                && informative_curvature
843            {
844                params.set_free_values(&free_vals_buf);
845                converged = true;
846                break;
847            }
848
849            // Reject step, restore parameters.
850            // y_current stays valid (parameters reverted to free_vals_buf snapshot).
851            params.set_free_values(&free_vals_buf);
852            lambda *= config.lambda_up;
853
854            // #108.4: If lambda is astronomically large, the optimizer is stuck
855            // in a region where no step improves chi2.  Break out rather than
856            // wasting iterations.
857            if lambda > LAMBDA_BREAKOUT {
858                converged = false;
859                break;
860            }
861        }
862    }
863
864    let reduced_chi2 = if dof > 0 { chi2 / dof as f64 } else { f64::NAN };
865
866    // Compute covariance matrix: (JᵀWJ)⁻¹ at the final parameters.
867    //
868    // This block requires an extra Jacobian evaluation + O(n_free³) matrix
869    // inversion.  When `compute_covariance` is false (e.g. per-pixel spatial
870    // mapping), we skip it entirely — the caller only needs densities and
871    // chi-squared, not uncertainties.
872    let (covariance, uncertainties) = if converged && config.compute_covariance {
873        let jacobian = compute_jacobian(
874            model,
875            params,
876            &y_current,
877            config.fd_step,
878            &mut all_vals_buf,
879            &mut free_idx_buf,
880        )?;
881        // Same active-mask row-skip semantics as the main assembly
882        // loop above — keeps the covariance free of NaN contributions
883        // from masked margin bins (#514).
884        let mut jtw_j = FlatMatrix::zeros(n_free, n_free);
885        for (i, &wi) in weights.iter().enumerate() {
886            if active_mask.is_some_and(|m| !m[i]) {
887                continue;
888            }
889            for j in 0..n_free {
890                let jij = jacobian.get(i, j);
891                for k in 0..n_free {
892                    *jtw_j.get_mut(j, k) += jij * wi * jacobian.get(i, k);
893                }
894            }
895        }
896
897        // #108.1: Scale covariance by reduced chi-squared.
898        //
899        // The raw (JᵀWJ)⁻¹ gives the covariance only when the model is a perfect
900        // description and the weights are exact.  Multiplying by χ²/ν accounts for
901        // misfit (model inadequacy or underestimated errors).  This is the standard
902        // statistical prescription (see e.g. Numerical Recipes §15.6).
903        //
904        // When dof == 0 (exactly determined system), reduced chi-squared is
905        // undefined (0/0).  We report NaN and skip covariance scaling entirely,
906        // returning None for covariance and uncertainties.
907        if dof > 0 {
908            if let Some(mut cov) = invert_matrix(&jtw_j) {
909                for elem in cov.data.iter_mut() {
910                    *elem *= reduced_chi2;
911                }
912                let unc: Vec<f64> = (0..n_free)
913                    .map(|i| {
914                        let diag = cov.get(i, i);
915                        if diag.is_finite() && diag > 0.0 {
916                            diag.sqrt()
917                        } else {
918                            f64::NAN
919                        }
920                    })
921                    .collect();
922                (Some(cov), Some(unc))
923            } else {
924                (None, None)
925            }
926        } else {
927            // dof == 0: covariance scaling is undefined; report None.
928            (None, None)
929        }
930    } else {
931        // Covariance computation skipped (compute_covariance == false).
932        (None, None)
933    };
934
935    Ok(LmResult {
936        chi_squared: chi2,
937        reduced_chi_squared: reduced_chi2,
938        iterations: iter,
939        converged,
940        params: params.all_values(),
941        covariance,
942        uncertainties,
943    })
944}
945
946#[cfg(test)]
947mod tests {
948    use super::*;
949    use crate::parameters::{FitParameter, ParameterSet};
950
951    /// Simple linear model: y = a*x + b
952    struct LinearModel {
953        x: Vec<f64>,
954    }
955
956    impl FitModel for LinearModel {
957        fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
958            let a = params[0];
959            let b = params[1];
960            Ok(self.x.iter().map(|&x| a * x + b).collect())
961        }
962    }
963
964    #[test]
965    fn test_fit_linear_exact() {
966        // Fit y = 2x + 3 with exact data (no noise)
967        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
968        let y_obs: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
969        let sigma = vec![1.0; 10];
970
971        let model = LinearModel { x };
972        let mut params = ParameterSet::new(vec![
973            FitParameter::unbounded("a", 1.0), // initial guess
974            FitParameter::unbounded("b", 1.0),
975        ]);
976
977        let result =
978            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
979
980        assert!(result.converged, "Fit did not converge");
981        assert!(
982            (result.params[0] - 2.0).abs() < 1e-4,
983            "a = {}, expected 2.0",
984            result.params[0]
985        );
986        assert!(
987            (result.params[1] - 3.0).abs() < 1e-4,
988            "b = {}, expected 3.0",
989            result.params[1]
990        );
991        assert!(result.chi_squared < 1e-6);
992    }
993
994    #[test]
995    fn test_converges_on_exact_flat_bottom_without_lambda_breakout() {
996        // Exact data with zero initial damping reaches the optimum in one
997        // Newton step. The next iteration sits on a flat χ² floor where the
998        // strict `trial_chi2 < chi2` check must not force a false non-
999        // convergence.
1000        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1001        let y_obs: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1002        let sigma = vec![1.0; 10];
1003
1004        let model = LinearModel { x };
1005        let mut params = ParameterSet::new(vec![
1006            FitParameter::unbounded("a", 0.0),
1007            FitParameter::unbounded("b", 0.0),
1008        ]);
1009        let config = LmConfig {
1010            lambda_init: 0.0,
1011            ..LmConfig::default()
1012        };
1013
1014        let result = levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &config).unwrap();
1015
1016        assert!(
1017            result.converged,
1018            "Fit should converge on an exact flat bottom"
1019        );
1020        assert!(result.chi_squared < 1e-20, "chi2 = {}", result.chi_squared);
1021        assert!(
1022            result.iterations < config.max_iter,
1023            "LM should stop by convergence, not by iteration exhaustion"
1024        );
1025    }
1026
1027    #[test]
1028    fn test_converges_on_nonzero_chi2_stationary_point() {
1029        // Noisy overdetermined data have a nonzero-chi2 optimum. With very
1030        // strict tolerances, the accepted step to the optimum does not trip
1031        // the accept-branch convergence checks, so the next iteration reaches
1032        // a reject-path stationary point with trial_chi2 == chi2.
1033        struct AffineModel {
1034            x: Vec<f64>,
1035        }
1036        impl FitModel for AffineModel {
1037            fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1038                let a = params[0];
1039                let b = params[1];
1040                Ok(self.x.iter().map(|&x| a * x + b).collect())
1041            }
1042        }
1043
1044        let model = AffineModel {
1045            x: vec![0.0, 1.0, 2.0, 3.0, 4.0],
1046        };
1047        let y_obs = vec![0.1, 0.9, 2.2, 2.8, 4.1];
1048        let sigma = vec![1.0; y_obs.len()];
1049        let mut params = ParameterSet::new(vec![
1050            FitParameter::unbounded("a", 0.0),
1051            FitParameter::unbounded("b", 0.0),
1052        ]);
1053        let config = LmConfig {
1054            max_iter: 200,
1055            tol_chi2: 1e-16,
1056            tol_param: 1e-16,
1057            ..LmConfig::default()
1058        };
1059
1060        let result = levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &config).unwrap();
1061
1062        assert!(
1063            result.converged,
1064            "stationary nonzero-chi2 optimum should converge instead of lambda breakout"
1065        );
1066        assert!(
1067            result.reduced_chi_squared.is_finite() && result.reduced_chi_squared > 0.0,
1068            "expected nonzero reduced chi2 at noisy optimum, got {}",
1069            result.reduced_chi_squared
1070        );
1071    }
1072
1073    #[test]
1074    fn test_fit_linear_with_fixed_param() {
1075        // Fit y = a*x + 3 with b fixed at 3.0
1076        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1077        let y_obs: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1078        let sigma = vec![1.0; 10];
1079
1080        let model = LinearModel { x };
1081        let mut params = ParameterSet::new(vec![
1082            FitParameter::unbounded("a", 1.0),
1083            FitParameter::fixed("b", 3.0),
1084        ]);
1085
1086        let result =
1087            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
1088
1089        assert!(result.converged);
1090        assert!(
1091            (result.params[0] - 2.0).abs() < 1e-6,
1092            "a = {}",
1093            result.params[0]
1094        );
1095        assert_eq!(result.params[1], 3.0); // fixed
1096    }
1097
1098    #[test]
1099    fn test_fit_quadratic() {
1100        // Fit y = a*x² + b*x + c to quadratic data
1101        struct QuadModel {
1102            x: Vec<f64>,
1103        }
1104        impl FitModel for QuadModel {
1105            fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1106                let (a, b, c) = (params[0], params[1], params[2]);
1107                Ok(self.x.iter().map(|&x| a * x * x + b * x + c).collect())
1108            }
1109        }
1110
1111        let x: Vec<f64> = (0..20).map(|i| i as f64 * 0.5).collect();
1112        let y_obs: Vec<f64> = x.iter().map(|&xi| 0.5 * xi * xi - 2.0 * xi + 1.0).collect();
1113        let sigma = vec![1.0; 20];
1114
1115        let model = QuadModel { x };
1116        let mut params = ParameterSet::new(vec![
1117            FitParameter::unbounded("a", 1.0),
1118            FitParameter::unbounded("b", 0.0),
1119            FitParameter::unbounded("c", 0.0),
1120        ]);
1121
1122        let result =
1123            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
1124
1125        assert!(result.converged);
1126        assert!(
1127            (result.params[0] - 0.5).abs() < 1e-5,
1128            "a = {}",
1129            result.params[0]
1130        );
1131        assert!(
1132            (result.params[1] - (-2.0)).abs() < 1e-5,
1133            "b = {}",
1134            result.params[1]
1135        );
1136        assert!(
1137            (result.params[2] - 1.0).abs() < 1e-5,
1138            "c = {}",
1139            result.params[2]
1140        );
1141    }
1142
1143    #[test]
1144    fn test_non_negative_constraint() {
1145        // Fit y = a*x with data that has negative slope,
1146        // but parameter a is constrained to be non-negative.
1147        // Should converge to a ≈ 0.
1148        struct SlopeModel {
1149            x: Vec<f64>,
1150        }
1151        impl FitModel for SlopeModel {
1152            fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1153                let a = params[0];
1154                Ok(self.x.iter().map(|&x| a * x).collect())
1155            }
1156        }
1157
1158        let x: Vec<f64> = (1..10).map(|i| i as f64).collect();
1159        let y_obs: Vec<f64> = x.iter().map(|&xi| -2.0 * xi).collect();
1160        let sigma = vec![1.0; 9];
1161
1162        let model = SlopeModel { x };
1163        let mut params = ParameterSet::new(vec![FitParameter::non_negative("a", 1.0)]);
1164
1165        let result =
1166            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
1167
1168        // Should be clamped at 0
1169        assert!(
1170            result.params[0] >= 0.0 && result.params[0] < 0.1,
1171            "a = {}, expected ~0",
1172            result.params[0]
1173        );
1174    }
1175
1176    #[test]
1177    fn test_uncertainty_estimation() {
1178        // Fit linear model; uncertainties should be reasonable
1179        let x: Vec<f64> = (0..100).map(|i| i as f64 * 0.1).collect();
1180        let y_obs: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1181        let sigma = vec![0.1; 100]; // Small uncertainty
1182
1183        let model = LinearModel { x };
1184        let mut params = ParameterSet::new(vec![
1185            FitParameter::unbounded("a", 1.0),
1186            FitParameter::unbounded("b", 1.0),
1187        ]);
1188
1189        let result =
1190            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
1191
1192        assert!(result.converged);
1193        assert!(result.uncertainties.is_some());
1194        let unc = result.uncertainties.unwrap();
1195        // Uncertainties should be positive and small
1196        assert!(unc[0] > 0.0 && unc[0] < 0.01, "σ_a = {}", unc[0]);
1197        assert!(unc[1] > 0.0 && unc[1] < 0.1, "σ_b = {}", unc[1]);
1198    }
1199
1200    #[test]
1201    fn test_solve_damped_system_identity() {
1202        // (I + λ·I)x = b → x = b/(1+λ)
1203        let a = FlatMatrix {
1204            data: vec![1.0, 0.0, 0.0, 1.0],
1205            nrows: 2,
1206            ncols: 2,
1207        };
1208        let b = vec![2.0, 4.0];
1209        let lambda = 1.0;
1210        let x = solve_damped_system(&a, &b, lambda).unwrap();
1211        assert!((x[0] - 1.0).abs() < 1e-10);
1212        assert!((x[1] - 2.0).abs() < 1e-10);
1213    }
1214
1215    #[test]
1216    fn test_invert_matrix_2x2() {
1217        let a = FlatMatrix {
1218            data: vec![4.0, 7.0, 2.0, 6.0],
1219            nrows: 2,
1220            ncols: 2,
1221        };
1222        let inv = invert_matrix(&a).unwrap();
1223        // A⁻¹ = 1/10 × [6 -7; -2 4]
1224        assert!((inv.get(0, 0) - 0.6).abs() < 1e-10);
1225        assert!((inv.get(0, 1) - (-0.7)).abs() < 1e-10);
1226        assert!((inv.get(1, 0) - (-0.2)).abs() < 1e-10);
1227        assert!((inv.get(1, 1) - 0.4).abs() < 1e-10);
1228    }
1229
1230    // ---- Edge-case tests for issue #125 ----
1231
1232    #[test]
1233    fn test_all_fixed_params_nan_model() {
1234        // #125.1: When all parameters are fixed and the model produces NaN,
1235        // the result must report converged=false (not converged=true with NaN chi2).
1236        struct NanModel;
1237        impl FitModel for NanModel {
1238            fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
1239                Ok(vec![f64::NAN; 5])
1240            }
1241        }
1242
1243        let y_obs = vec![1.0; 5];
1244        let sigma = vec![1.0; 5];
1245        let mut params = ParameterSet::new(vec![FitParameter::fixed("a", 1.0)]);
1246
1247        let result =
1248            levenberg_marquardt(&NanModel, &y_obs, &sigma, &mut params, &LmConfig::default())
1249                .unwrap();
1250
1251        assert!(!result.converged, "All-fixed NaN model should not converge");
1252        assert!(result.chi_squared.is_nan(), "chi2 should be NaN");
1253        assert_eq!(result.iterations, 0);
1254    }
1255
1256    #[test]
1257    fn test_underdetermined_system() {
1258        // #125.6: More free parameters than data points → underdetermined.
1259        // Should return converged=false immediately.
1260        let y_obs = vec![1.0, 2.0]; // 2 data points
1261        let sigma = vec![1.0, 1.0];
1262
1263        // 2 free params for 2 data points is exactly determined (ok),
1264        // but 3 free params for 2 data points is underdetermined.
1265        struct ThreeParamModel;
1266        impl FitModel for ThreeParamModel {
1267            fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1268                Ok(vec![params[0] + params[1] + params[2]; 2])
1269            }
1270        }
1271
1272        let mut params = ParameterSet::new(vec![
1273            FitParameter::unbounded("a", 1.0),
1274            FitParameter::unbounded("b", 1.0),
1275            FitParameter::unbounded("c", 1.0),
1276        ]);
1277
1278        let result = levenberg_marquardt(
1279            &ThreeParamModel,
1280            &y_obs,
1281            &sigma,
1282            &mut params,
1283            &LmConfig::default(),
1284        )
1285        .unwrap();
1286
1287        assert!(
1288            !result.converged,
1289            "Underdetermined system should not converge"
1290        );
1291        assert!(result.chi_squared.is_nan());
1292        assert_eq!(result.iterations, 0);
1293    }
1294
1295    #[test]
1296    fn test_exactly_determined_dof_zero() {
1297        // #125.6: n_data == n_free → dof=0, exactly determined.
1298        // Should still converge but reduced_chi_squared is NaN (0/0).
1299        let y_obs = vec![5.0, 11.0]; // y = 2x + 3 at x=1,4
1300        let sigma = vec![1.0, 1.0];
1301
1302        let model = LinearModel { x: vec![1.0, 4.0] };
1303        let mut params = ParameterSet::new(vec![
1304            FitParameter::unbounded("a", 1.0),
1305            FitParameter::unbounded("b", 1.0),
1306        ]);
1307
1308        let result =
1309            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
1310
1311        assert!(
1312            result.converged,
1313            "Exactly-determined system should converge"
1314        );
1315        assert!(
1316            result.chi_squared < 1e-6,
1317            "chi2 should be ~0, got {}",
1318            result.chi_squared
1319        );
1320        // dof=0 → reduced chi2 is NaN
1321        assert!(
1322            result.reduced_chi_squared.is_nan(),
1323            "reduced_chi2 should be NaN for dof=0, got {}",
1324            result.reduced_chi_squared
1325        );
1326        // No covariance when dof=0
1327        assert!(result.covariance.is_none());
1328        assert!(result.uncertainties.is_none());
1329    }
1330
1331    #[test]
1332    fn test_lambda_breakout() {
1333        // #125.6: A model that never improves should trigger lambda breakout.
1334        struct ConstantModel;
1335        impl FitModel for ConstantModel {
1336            fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
1337                // Returns constant output regardless of parameters,
1338                // so the Jacobian is zero and no step can improve chi2.
1339                Ok(vec![42.0; 5])
1340            }
1341        }
1342
1343        let y_obs = vec![1.0, 2.0, 3.0, 4.0, 5.0];
1344        let sigma = vec![1.0; 5];
1345        let mut params = ParameterSet::new(vec![FitParameter::unbounded("a", 1.0)]);
1346
1347        let config = LmConfig {
1348            max_iter: 1000,
1349            ..LmConfig::default()
1350        };
1351
1352        let result =
1353            levenberg_marquardt(&ConstantModel, &y_obs, &sigma, &mut params, &config).unwrap();
1354
1355        assert!(
1356            !result.converged,
1357            "Flat model should not converge (lambda breakout)"
1358        );
1359        assert!(
1360            result.covariance.is_none(),
1361            "unconverged fit should not report covariance"
1362        );
1363        assert!(
1364            result.uncertainties.is_none(),
1365            "unconverged fit should not report uncertainties"
1366        );
1367    }
1368
1369    #[test]
1370    fn test_nan_model_during_iteration() {
1371        // #125.6: Model that produces NaN for certain parameter values.
1372        // The optimizer should treat NaN steps as bad and try smaller steps.
1373        struct NanAtLargeModel {
1374            x: Vec<f64>,
1375        }
1376        impl FitModel for NanAtLargeModel {
1377            fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1378                let a = params[0];
1379                Ok(self
1380                    .x
1381                    .iter()
1382                    .map(|&x| if a > 5.0 { f64::NAN } else { a * x + 1.0 })
1383                    .collect())
1384            }
1385        }
1386
1387        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1388        let y_obs: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 1.0).collect();
1389        let sigma = vec![1.0; 10];
1390
1391        let model = NanAtLargeModel { x };
1392        let mut params = ParameterSet::new(vec![FitParameter::unbounded("a", 3.0)]);
1393
1394        let result =
1395            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
1396
1397        // Should converge to a≈2 while avoiding the NaN region a>5.
1398        assert!(result.converged, "Should converge avoiding NaN region");
1399        assert!(
1400            (result.params[0] - 2.0).abs() < 0.1,
1401            "a = {}, expected ~2.0",
1402            result.params[0]
1403        );
1404    }
1405
1406    #[test]
1407    fn test_err_model_during_trial_step() {
1408        // Model that returns Err for large parameter values.
1409        // The optimizer should treat Err trial steps as bad steps (increase λ)
1410        // and converge without panicking.
1411        struct ErrAtLargeModel {
1412            x: Vec<f64>,
1413        }
1414        impl FitModel for ErrAtLargeModel {
1415            fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1416                let a = params[0];
1417                if a > 5.0 {
1418                    return Err(FittingError::EvaluationFailed(
1419                        "parameter out of valid range".into(),
1420                    ));
1421                }
1422                Ok(self.x.iter().map(|&x| a * x + 1.0).collect())
1423            }
1424        }
1425
1426        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1427        let y_obs: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 1.0).collect();
1428        let sigma = vec![1.0; 10];
1429
1430        let model = ErrAtLargeModel { x };
1431        let mut params = ParameterSet::new(vec![FitParameter::unbounded("a", 3.0)]);
1432
1433        let result =
1434            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
1435
1436        // Should converge to a≈2 while avoiding the Err region a>5.
1437        assert!(result.converged, "Should converge avoiding Err region");
1438        assert!(
1439            (result.params[0] - 2.0).abs() < 0.1,
1440            "a = {}, expected ~2.0",
1441            result.params[0]
1442        );
1443    }
1444
1445    #[test]
1446    fn test_fit_linear_no_covariance() {
1447        // When compute_covariance is false, the fit should still converge and
1448        // produce correct parameters, but covariance and uncertainties are None.
1449        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1450        let y_obs: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1451        let sigma = vec![1.0; 10];
1452
1453        let model = LinearModel { x };
1454        let mut params = ParameterSet::new(vec![
1455            FitParameter::unbounded("a", 1.0),
1456            FitParameter::unbounded("b", 1.0),
1457        ]);
1458
1459        let config = LmConfig {
1460            compute_covariance: false,
1461            ..LmConfig::default()
1462        };
1463
1464        let result = levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &config).unwrap();
1465
1466        assert!(result.converged, "Fit did not converge");
1467        assert!(
1468            (result.params[0] - 2.0).abs() < 1e-4,
1469            "a = {}, expected 2.0",
1470            result.params[0]
1471        );
1472        assert!(
1473            (result.params[1] - 3.0).abs() < 1e-4,
1474            "b = {}, expected 3.0",
1475            result.params[1]
1476        );
1477        assert!(result.chi_squared < 1e-6);
1478        assert!(
1479            result.covariance.is_none(),
1480            "covariance should be None when compute_covariance=false"
1481        );
1482        assert!(
1483            result.uncertainties.is_none(),
1484            "uncertainties should be None when compute_covariance=false"
1485        );
1486    }
1487
1488    #[test]
1489    fn test_zero_negative_sigma_clamping() {
1490        // #125.6: Zero and negative sigma should be clamped to huge sigma (tiny weight),
1491        // not cause NaN/panic.
1492        let y_obs = vec![1.0, 2.0, 3.0, 4.0, 5.0];
1493        let sigma = vec![0.0, -1.0, f64::NAN, f64::INFINITY, 1.0];
1494
1495        let model = LinearModel {
1496            x: vec![0.0, 1.0, 2.0, 3.0, 4.0],
1497        };
1498        let mut params = ParameterSet::new(vec![
1499            FitParameter::unbounded("a", 1.0),
1500            FitParameter::unbounded("b", 0.0),
1501        ]);
1502
1503        // Should not panic and should produce a finite result.
1504        let result =
1505            levenberg_marquardt(&model, &y_obs, &sigma, &mut params, &LmConfig::default()).unwrap();
1506
1507        assert!(
1508            result.chi_squared.is_finite(),
1509            "chi2 should be finite despite bad sigma, got {}",
1510            result.chi_squared
1511        );
1512        assert!(
1513            result.converged,
1514            "Fit should converge despite bad sigma values"
1515        );
1516        // The only valid data point with sigma=1.0 is (x=4, y=5).
1517        // The fitted line y = a*x + b should pass near that point.
1518        let y_at_4 = result.params[0] * 4.0 + result.params[1];
1519        assert!(
1520            (y_at_4 - 5.0).abs() < 1.0,
1521            "Fitted line should pass near (4, 5): a={}, b={}, y(4)={}",
1522            result.params[0],
1523            result.params[1],
1524            y_at_4,
1525        );
1526    }
1527
1528    // ------------------------------------------------------------------
1529    // Active-bin mask (SAMMY EMIN/EMAX-equivalent fit-energy-range, #514).
1530    // ------------------------------------------------------------------
1531
1532    /// LM with an `active_mask` that restricts to a subset of bins must
1533    /// produce a fit equivalent to running on that subset directly —
1534    /// validates that masking propagates through both the residual /
1535    /// χ² accumulation AND the dof / reduced-χ² book-keeping.
1536    #[test]
1537    fn test_lm_active_mask_subset_equivalence() {
1538        // Linear fit y = 2x + 3 on x = 0..10.  Mask = bins 0..5 only.
1539        let x_full: Vec<f64> = (0..10).map(|i| i as f64).collect();
1540        let y_full: Vec<f64> = x_full.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1541        let sigma_full = vec![1.0; 10];
1542        let mask = vec![
1543            true, true, true, true, true, false, false, false, false, false,
1544        ];
1545
1546        let model = LinearModel { x: x_full.clone() };
1547        let mut params = ParameterSet::new(vec![
1548            FitParameter::unbounded("a", 1.0),
1549            FitParameter::unbounded("b", 1.0),
1550        ]);
1551        let result_masked = levenberg_marquardt_with_mask(
1552            &model,
1553            &y_full,
1554            &sigma_full,
1555            &mut params,
1556            &LmConfig::default(),
1557            Some(&mask),
1558        )
1559        .unwrap();
1560        assert!(result_masked.converged);
1561        assert!((result_masked.params[0] - 2.0).abs() < 1e-6);
1562        assert!((result_masked.params[1] - 3.0).abs() < 1e-6);
1563
1564        // Reference: same fit on the subset directly (no mask).
1565        let x_sub = x_full[..5].to_vec();
1566        let y_sub = y_full[..5].to_vec();
1567        let sigma_sub = vec![1.0; 5];
1568        let model_sub = LinearModel { x: x_sub };
1569        let mut params_sub = ParameterSet::new(vec![
1570            FitParameter::unbounded("a", 1.0),
1571            FitParameter::unbounded("b", 1.0),
1572        ]);
1573        let result_sub = levenberg_marquardt(
1574            &model_sub,
1575            &y_sub,
1576            &sigma_sub,
1577            &mut params_sub,
1578            &LmConfig::default(),
1579        )
1580        .unwrap();
1581
1582        // Recovered parameters and reduced-χ² (both fits noiseless,
1583        // so reduced-χ² ≈ 0 in either case) must match.
1584        for j in 0..2 {
1585            assert!(
1586                (result_masked.params[j] - result_sub.params[j]).abs() < 1e-6,
1587                "param {j}: masked={} subset={}",
1588                result_masked.params[j],
1589                result_sub.params[j]
1590            );
1591        }
1592        assert!(
1593            (result_masked.reduced_chi_squared - result_sub.reduced_chi_squared).abs() < 1e-6,
1594            "reduced-χ²: masked={} subset={}",
1595            result_masked.reduced_chi_squared,
1596            result_sub.reduced_chi_squared
1597        );
1598    }
1599
1600    /// Without the mask, residuals from the out-of-range half of the
1601    /// grid would dominate χ² and pull the fit way off — the masked
1602    /// fit must NOT be biased by them.
1603    #[test]
1604    fn test_lm_active_mask_excludes_out_of_range_residuals() {
1605        // y = 2x + 3 inside the masked range; corrupt outliers outside.
1606        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1607        let mut y: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1608        for yi in y.iter_mut().skip(5) {
1609            *yi += 1000.0; // huge corruption outside the active mask
1610        }
1611        let sigma = vec![1.0; 10];
1612        let mask = vec![true; 5]
1613            .into_iter()
1614            .chain(std::iter::repeat_n(false, 5))
1615            .collect::<Vec<bool>>();
1616
1617        let model = LinearModel { x };
1618        let mut params = ParameterSet::new(vec![
1619            FitParameter::unbounded("a", 1.0),
1620            FitParameter::unbounded("b", 1.0),
1621        ]);
1622        let result = levenberg_marquardt_with_mask(
1623            &model,
1624            &y,
1625            &sigma,
1626            &mut params,
1627            &LmConfig::default(),
1628            Some(&mask),
1629        )
1630        .unwrap();
1631        assert!(result.converged);
1632        assert!(
1633            (result.params[0] - 2.0).abs() < 1e-6,
1634            "slope should be 2.0 (corrupt outliers must be masked out), got {}",
1635            result.params[0]
1636        );
1637        assert!(
1638            (result.params[1] - 3.0).abs() < 1e-6,
1639            "intercept should be 3.0 (corrupt outliers must be masked out), got {}",
1640            result.params[1]
1641        );
1642    }
1643
1644    /// Out-of-mask `y_obs` containing `NaN` must not poison χ² /
1645    /// JᵀWJ / JᵀWr.  Without explicit row-skip in the accumulator
1646    /// loops, `0.0 (weight) * NaN (residual) = NaN` would propagate
1647    /// through the fit despite the masked weight being zero.
1648    /// Regression for Round-2 review fix #4.
1649    #[test]
1650    fn test_lm_active_mask_tolerates_nan_outside_range() {
1651        // y = 2x + 3 inside the active mask; NaN outside.  A naive
1652        // zero-weight implementation would return NaN χ² and fail to
1653        // converge; the row-skip path should fit cleanly.
1654        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1655        let mut y: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1656        for yi in y.iter_mut().skip(5) {
1657            *yi = f64::NAN;
1658        }
1659        let sigma = vec![1.0; 10];
1660        let mask = vec![true; 5]
1661            .into_iter()
1662            .chain(std::iter::repeat_n(false, 5))
1663            .collect::<Vec<bool>>();
1664
1665        let model = LinearModel { x };
1666        let mut params = ParameterSet::new(vec![
1667            FitParameter::unbounded("a", 1.0),
1668            FitParameter::unbounded("b", 1.0),
1669        ]);
1670        let result = levenberg_marquardt_with_mask(
1671            &model,
1672            &y,
1673            &sigma,
1674            &mut params,
1675            &LmConfig::default(),
1676            Some(&mask),
1677        )
1678        .unwrap();
1679
1680        assert!(
1681            result.converged,
1682            "fit should converge despite NaN outside the active mask"
1683        );
1684        assert!(
1685            result.chi_squared.is_finite(),
1686            "χ² should be finite (NaN poisoning prevented), got {}",
1687            result.chi_squared
1688        );
1689        assert!(
1690            (result.params[0] - 2.0).abs() < 1e-6,
1691            "slope should be 2.0, got {}",
1692            result.params[0]
1693        );
1694        assert!(
1695            (result.params[1] - 3.0).abs() < 1e-6,
1696            "intercept should be 3.0, got {}",
1697            result.params[1]
1698        );
1699    }
1700
1701    /// Regression for Round-2-pass-2 review: the global finite checks
1702    /// on `y_model` / `y_trial` must skip masked bins.  A model that
1703    /// returns NaN only at masked margin bins should not abort the fit
1704    /// or get its trial step rejected.
1705    #[test]
1706    fn test_lm_active_mask_tolerates_model_nan_outside_range() {
1707        // Linear model that injects NaN into masked bins.  The fit
1708        // should still converge cleanly — the contract is that masked
1709        // rows are completely irrelevant to the fit.
1710        struct LinearModelWithMaskedNaN {
1711            x: Vec<f64>,
1712            mask: Vec<bool>,
1713        }
1714        impl FitModel for LinearModelWithMaskedNaN {
1715            fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1716                let a = params[0];
1717                let b = params[1];
1718                Ok(self
1719                    .x
1720                    .iter()
1721                    .enumerate()
1722                    .map(|(i, &x)| if self.mask[i] { a * x + b } else { f64::NAN })
1723                    .collect())
1724            }
1725        }
1726
1727        let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
1728        let y: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1729        let sigma = vec![1.0; 10];
1730        let mask: Vec<bool> = vec![true; 5]
1731            .into_iter()
1732            .chain(std::iter::repeat_n(false, 5))
1733            .collect();
1734
1735        let model = LinearModelWithMaskedNaN {
1736            x: x.clone(),
1737            mask: mask.clone(),
1738        };
1739        let mut params = ParameterSet::new(vec![
1740            FitParameter::unbounded("a", 1.0),
1741            FitParameter::unbounded("b", 1.0),
1742        ]);
1743        let result = levenberg_marquardt_with_mask(
1744            &model,
1745            &y,
1746            &sigma,
1747            &mut params,
1748            &LmConfig::default(),
1749            Some(&mask),
1750        )
1751        .unwrap();
1752
1753        assert!(
1754            result.converged,
1755            "fit should converge despite model NaN at masked bins"
1756        );
1757        assert!(result.chi_squared.is_finite());
1758        assert!((result.params[0] - 2.0).abs() < 1e-6);
1759        assert!((result.params[1] - 3.0).abs() < 1e-6);
1760    }
1761
1762    /// A zero-active mask (the user's range misses the grid entirely)
1763    /// must return non-converged regardless of whether parameters are
1764    /// free or fixed.  Pre-fix the all-fixed fast-return path would
1765    /// report `converged: true, chi_squared: 0` from a sum over zero
1766    /// rows — a deceptive "success" with no data behind it.
1767    #[test]
1768    fn test_lm_active_mask_all_false_returns_non_converged() {
1769        let x: Vec<f64> = (0..5).map(|i| i as f64).collect();
1770        let y: Vec<f64> = x.iter().map(|&xi| 2.0 * xi + 3.0).collect();
1771        let sigma = vec![1.0; 5];
1772        let mask = vec![false; 5];
1773        let model = LinearModel { x: x.clone() };
1774
1775        // (a) All parameters free.
1776        let mut params_free = ParameterSet::new(vec![
1777            FitParameter::unbounded("a", 1.0),
1778            FitParameter::unbounded("b", 1.0),
1779        ]);
1780        let r_free = levenberg_marquardt_with_mask(
1781            &model,
1782            &y,
1783            &sigma,
1784            &mut params_free,
1785            &LmConfig::default(),
1786            Some(&mask),
1787        )
1788        .unwrap();
1789        assert!(
1790            !r_free.converged,
1791            "n_free > 0 + zero-active mask must NOT report converged"
1792        );
1793        assert!(r_free.chi_squared.is_nan());
1794        assert!(r_free.reduced_chi_squared.is_nan());
1795
1796        // (b) All parameters fixed → n_free == 0 (the path #517 R3
1797        //     specifically failed before the n_active==0 early-return).
1798        let mut params_fixed = ParameterSet::new(vec![
1799            FitParameter::fixed("a", 1.0),
1800            FitParameter::fixed("b", 0.0),
1801        ]);
1802        let r_fixed = levenberg_marquardt_with_mask(
1803            &model,
1804            &y,
1805            &sigma,
1806            &mut params_fixed,
1807            &LmConfig::default(),
1808            Some(&mask),
1809        )
1810        .unwrap();
1811        assert!(
1812            !r_fixed.converged,
1813            "n_free == 0 + zero-active mask must NOT report converged \
1814             (sum over zero rows would be 0, masquerading as a perfect fit)"
1815        );
1816        assert!(r_fixed.chi_squared.is_nan());
1817        assert!(r_fixed.reduced_chi_squared.is_nan());
1818    }
1819
1820    // ==================================================================
1821    // NaN-in-Jacobian during FD probes.
1822    //
1823    // The main LM loop guards the trial step at the `model.evaluate`
1824    // site via `trial_has_active_nonfinite`.  The silent surface is the
1825    // post-convergence covariance Jacobian: it calls `compute_jacobian`
1826    // directly, which has no finiteness check on the per-column FD
1827    // probe.  A NaN at any active row of the perturbed model output
1828    // gets divided by `actual_step` and turned into NaN entries in the
1829    // Jacobian — these poison `JᵀWJ` and the inverse covariance, which
1830    // is reported as the "uncertainty" of the fit.
1831    //
1832    // `compute_jacobian`'s FD path zeros per-cell entries whose probe
1833    // returned a non-finite value, mirroring the existing `Err` branch.
1834    // ==================================================================
1835
1836    /// `compute_jacobian`'s FD path zeroes per-cell entries whose
1837    /// perturbed model output is non-finite, rather than baking NaN
1838    /// into the Jacobian (which then poisons the post-convergence
1839    /// covariance computation).  Per-cell (not whole-column) so that a
1840    /// NaN at a masked / inactive row leaves the rest of the column
1841    /// usable — see `test_lm_active_mask_tolerates_model_nan_outside_range`
1842    /// for the masked-NaN contract.
1843    #[test]
1844    fn test_compute_jacobian_skips_nan_perturbed_column() {
1845        use crate::parameters::FitParameter;
1846        // Two-parameter model.  Param 0 is well-behaved (returns
1847        // `θ_0` constant).  Param 1 produces a NaN row at the +FD
1848        // probe — exactly the "NaN-in-Jacobian during FD probes"
1849        // signature from issue #552.
1850        struct NanInColumn1 {
1851            // Base value of param 1 — used to detect the perturbation.
1852            p1_base: f64,
1853        }
1854        impl FitModel for NanInColumn1 {
1855            fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
1856                // The +FD probe perturbs params[1] by ε; everything else
1857                // keeps params[1] == self.p1_base.
1858                let nan_row = (params[1] - self.p1_base).abs() > 1e-12;
1859                let v = if nan_row { f64::NAN } else { params[0] };
1860                Ok(vec![v; 4])
1861            }
1862            // No analytical_jacobian -> FD fallback drives compute_jacobian.
1863        }
1864        let model = NanInColumn1 { p1_base: 0.3 };
1865        let mut params = ParameterSet::new(vec![
1866            FitParameter::unbounded("p0", 0.5),
1867            FitParameter::unbounded("p1", 0.3),
1868        ]);
1869        let y_current = vec![0.5; 4];
1870        let mut all_vals_buf: Vec<f64> = Vec::new();
1871        let mut free_idx_buf: Vec<usize> = Vec::new();
1872        let jac = compute_jacobian(
1873            &model,
1874            &mut params,
1875            &y_current,
1876            1e-6,
1877            &mut all_vals_buf,
1878            &mut free_idx_buf,
1879        )
1880        .unwrap();
1881        // Every entry must be finite — column 1 must have been skipped.
1882        for (i, v) in jac.data.iter().enumerate() {
1883            assert!(
1884                v.is_finite(),
1885                "compute_jacobian produced non-finite entry at index {i}: {v}"
1886            );
1887        }
1888        // Column 1 was skipped -> all zeros.  Column 0 has the normal
1889        // ∂T/∂θ_0 = 1 column from the well-behaved param.
1890        for i in 0..jac.nrows {
1891            assert_eq!(
1892                jac.get(i, 1),
1893                0.0,
1894                "column 1 (NaN probe) should be zeroed, row {i}"
1895            );
1896        }
1897    }
1898}