Skip to main content

nereids_physics/
transmission.rs

1//! Transmission forward model via the Beer-Lambert law.
2//!
3//! Computes theoretical neutron transmission spectra from resonance parameters,
4//! applying cross-section calculation, Doppler broadening, resolution broadening,
5//! and the Beer-Lambert attenuation law.
6//!
7//! ## Beer-Lambert Law
8//!
9//! For a single isotope:
10//!   T(E) = exp(-n·d·σ(E))
11//!
12//! For multiple isotopes:
13//!   T(E) = exp(-Σᵢ nᵢ·dᵢ·σᵢ(E))
14//!
15//! where n is number density (atoms/cm³), d is thickness (cm),
16//! and σ(E) is the total cross-section in barns (1 barn = 10⁻²⁴ cm²).
17//!
18//! In practice, the product n·d is expressed as "areal density" in
19//! atoms/barn, so T(E) = exp(-thickness × σ(E)) with thickness in atoms/barn.
20//!
21//! ## SAMMY Reference
22//! - `cro/` and `xxx/` modules — cross-section to transmission conversion
23//! - Manual Section II (cross-section theory); transmission experiments
24//!   Section III.E.1
25
26use std::fmt;
27use std::sync::atomic::{AtomicBool, Ordering};
28
29use rayon::prelude::*;
30
31use nereids_endf::resonance::ResonanceData;
32
33use crate::doppler::{self, DopplerParams, DopplerParamsError};
34use crate::reich_moore;
35use crate::resolution::{self, ResolutionError, ResolutionFunction};
36
37/// Build the auxiliary extended grid for resolution broadening.
38///
39/// Shared helper that extracts Gaussian resolution params and resonance info
40/// to build the extended grid with boundary extension + adaptive intermediate
41/// points.  Returns `None` if no extension is needed (no resolution, or grid
42/// unchanged).
43///
44/// Intermediate points are inserted only when the resolution broadening at
45/// the grid midpoint uses the PW-linear Gaussian path (exp tail negligible
46/// or absent).  For genuine combined-kernel cases, intermediates create
47/// non-uniform spacing transitions that degrade the Xcoef quadrature.
48fn build_aux_grid(
49    energies: &[f64],
50    instrument: Option<&InstrumentParams>,
51    resonance_data: &[&ResonanceData],
52) -> Option<(Vec<f64>, Vec<usize>)> {
53    instrument.and_then(|inst| {
54        if let ResolutionFunction::Gaussian(ref params) = inst.resolution {
55            // P-9: Check the Gaussian-to-exp-tail ratio at MULTIPLE energies
56            // to decide whether intermediates help or hurt.  The ratio C =
57            // W_g/(2·W_e) determines which broadening path is used per-energy
58            // in resolution_broaden_presorted.  When C > 2.5 the PW-linear
59            // Gaussian path is used, which benefits from intermediates.
60            //
61            // Previously checked at a single midpoint, which could make the
62            // wrong decision if the ratio crosses 2.5 within the energy range.
63            // Now checks at 5 points (lo, 25%, mid, 75%, hi) and uses
64            // intermediates if a MAJORITY of points have C > 2.5.
65            let use_intermediates = if energies.len() >= 2 {
66                let n = energies.len();
67                let check_indices = [0, n / 4, n / 2, 3 * n / 4, n - 1];
68                let n_pw_linear = check_indices
69                    .iter()
70                    .filter(|&&i| {
71                        let e = energies[i];
72                        let wg = params.gaussian_width(e);
73                        let we = params.exp_width(e);
74                        we < 1e-60 || wg / (2.0 * we) > 2.5
75                    })
76                    .count();
77                // Majority rule: use intermediates if ≥3 of 5 points qualify.
78                n_pw_linear >= 3
79            } else {
80                true
81            };
82
83            // Extract (energy_eV, gd_eV) pairs for fine-structure densification.
84            // gd = total resonance width, used by Fspken to identify regions
85            // needing denser grid points around narrow resonances.
86            // SAMMY Ref: dat/mdat4.f90 Fspken lines 243-284
87            let resonances = extract_resonance_widths(resonance_data);
88
89            let (ext_e, di) = if use_intermediates {
90                crate::auxiliary_grid::build_extended_grid(energies, Some(params), &resonances)
91            } else {
92                crate::auxiliary_grid::build_extended_grid_boundary_only(energies, Some(params))
93            };
94            if ext_e.len() > energies.len() {
95                Some((ext_e, di))
96            } else {
97                None
98            }
99        } else {
100            None
101        }
102    })
103}
104
105/// Extract (energy_eV, gd_eV) pairs from resonance data for fine-structure
106/// grid densification.
107///
108/// For LRF=1/2/3 (BW and Reich-Moore): `gd = |Γn| + |Γγ| + |Γf1| + |Γf2|`
109/// For LRF=7 (R-Matrix Limited): `gd = |Γγ| + Σ|γ_i|²` (approximate)
110///
111/// SAMMY Ref: dat/mdat4.f90 Fspken — uses total width to define the region
112/// [E_res − gd, E_res + gd] for fine-structure point insertion.
113fn extract_resonance_widths(resonance_data: &[&ResonanceData]) -> Vec<(f64, f64)> {
114    let mut pairs = Vec::new();
115    for rd in resonance_data {
116        for range in &rd.ranges {
117            if !range.resolved {
118                continue;
119            }
120            // LRF=1/2/3: resonances grouped by L
121            for lg in &range.l_groups {
122                for res in &lg.resonances {
123                    let gd = res.gn.abs() + res.gg.abs() + res.gfa.abs() + res.gfb.abs();
124                    if gd > 0.0 {
125                        pairs.push((res.energy, gd));
126                    }
127                }
128            }
129            // LRF=7: resonances in spin groups
130            if let Some(ref rml) = range.rml {
131                for sg in &rml.spin_groups {
132                    for res in &sg.resonances {
133                        let mut gd = res.gamma_gamma.abs();
134                        for &w in &res.widths {
135                            gd += w * w; // γ² approximates Γ
136                        }
137                        if gd > 0.0 {
138                            pairs.push((res.energy, gd));
139                        }
140                    }
141                }
142            }
143        }
144    }
145    pairs
146}
147
148/// Distinct resonance center energies (eV) across every isotope and formalism,
149/// sorted ascending.
150///
151/// Used by the energy-scale calibration peak-matching seed (issue #608): the
152/// (t0, L_scale) calibration is recovered by matching measured transmission-dip
153/// positions to these known resonance energies via the linear TOF relation
154/// `tof_measured = t0 + L_scale · tof_nominal`, which seeds the LM into the
155/// global-minimum basin of the (sharply non-convex, post-#608) calibration χ²
156/// surface — a basin too thin for a cold start or grid scan to find reliably.
157///
158/// **Deduplicated** (Copilot PR #609): two resonances at the same center energy
159/// (e.g. across isotopes in a grouped fit) are a single POSITION for dip
160/// matching — multiplicity is irrelevant — and a duplicate would drive the
161/// minimum inter-resonance spacing, and hence the seed's `match_tol`, to 0,
162/// rejecting every dip and silently dropping the seed to a cold start.
163pub fn resonance_center_energies(resonance_data: &[&ResonanceData]) -> Vec<f64> {
164    let mut e: Vec<f64> = extract_resonance_widths(resonance_data)
165        .into_iter()
166        .map(|(energy, _gd)| energy)
167        .collect();
168    e.sort_by(f64::total_cmp);
169    // Sorted ⇒ exact duplicates are adjacent; `dedup` removes them.
170    e.dedup();
171    e
172}
173
174/// Build the unbroadened cross-section vector on the extended/auxiliary grid
175/// from cached data-grid values.
176///
177/// Data-grid positions are copied verbatim from `xs_raw` (the cached
178/// unbroadened σ for this isotope); the auxiliary-only positions (boundary
179/// extension + fine-structure points) are evaluated fresh.  `is_data_point`
180/// marks which extended-grid indices are data-grid points.
181///
182/// This is the cheap reuse of cached XS that the base-XS family relies on: only
183/// the few hundred auxiliary points are recomputed, not the full grid.
184fn build_extended_xs_from_base(
185    ext_energies: &[f64],
186    data_indices: &[usize],
187    is_data_point: &[bool],
188    xs_raw: &[f64],
189    rd: &ResonanceData,
190) -> Vec<f64> {
191    let mut xs_ext = vec![0.0f64; ext_energies.len()];
192    for (data_i, &ext_i) in data_indices.iter().enumerate() {
193        xs_ext[ext_i] = xs_raw[data_i];
194    }
195    for (j, &e) in ext_energies.iter().enumerate() {
196        if !is_data_point[j] {
197            xs_ext[j] = reich_moore::cross_sections_at_energy(rd, e).total;
198        }
199    }
200    xs_ext
201}
202
203/// Validate `base_xs` shape against `resonance_data` and `energies`.
204///
205/// Shared by the base-XS broadening family so the same `InputMismatch`
206/// diagnostics are produced from one place.
207fn validate_base_xs(
208    energies: &[f64],
209    base_xs: &[Vec<f64>],
210    resonance_data: &[ResonanceData],
211) -> Result<(), TransmissionError> {
212    if base_xs.len() != resonance_data.len() {
213        return Err(TransmissionError::InputMismatch(format!(
214            "base_xs has {} isotopes but resonance_data has {}",
215            base_xs.len(),
216            resonance_data.len(),
217        )));
218    }
219    for (i, row) in base_xs.iter().enumerate() {
220        if row.len() != energies.len() {
221            return Err(TransmissionError::InputMismatch(format!(
222                "base_xs[{i}] has {} energies but expected {}",
223                row.len(),
224                energies.len(),
225            )));
226        }
227    }
228    Ok(())
229}
230
231/// Build a bool mask identifying data-grid positions within the extended grid.
232///
233/// `None` when no auxiliary grid was built (the working grid is the data grid).
234fn data_point_mask(ext_grid: Option<&(Vec<f64>, Vec<usize>)>) -> Option<Vec<bool>> {
235    ext_grid.map(|(ext_e, di)| {
236        let mut mask = vec![false; ext_e.len()];
237        for &idx in di {
238            mask[idx] = true;
239        }
240        mask
241    })
242}
243
244/// Resolve the working grid (extended grid when available, else the data grid)
245/// and the [`WorkingGridLayout`] mapping data points back into it.
246fn working_grid_layout<'a>(
247    energies: &'a [f64],
248    ext_grid: Option<&'a (Vec<f64>, Vec<usize>)>,
249) -> (&'a [f64], WorkingGridLayout) {
250    match ext_grid {
251        Some((ext_e, di)) => (
252            ext_e.as_slice(),
253            WorkingGridLayout {
254                energies: ext_e.clone(),
255                data_indices: di.clone(),
256            },
257        ),
258        None => (
259            energies,
260            WorkingGridLayout {
261                energies: energies.to_vec(),
262                data_indices: (0..energies.len()).collect(),
263            },
264        ),
265    }
266}
267
268/// Compute the working-grid layout for `(energies, instrument, resonance_data)`.
269///
270/// Returns the auxiliary extended grid (boundary extension + resonance
271/// fine-structure) when a Gaussian resolution function is active, otherwise the
272/// data grid with identity indices.  Fitting models use this to apply
273/// resolution on the same working grid [`forward_model`] uses and extract the
274/// data points last (issue #608).
275///
276/// `resonance_data` may be empty (e.g. the energy-scale model has no resonance
277/// data of its own); the auxiliary grid then carries boundary extension only,
278/// with no resonance fine-structure densification.
279///
280/// # Errors
281/// * [`TransmissionError::Resolution`] — if `instrument` is `Some` and
282///   `energies` is not sorted ascending.
283pub fn resolution_working_grid(
284    energies: &[f64],
285    instrument: Option<&InstrumentParams>,
286    resonance_data: &[&ResonanceData],
287) -> Result<WorkingGridLayout, TransmissionError> {
288    if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
289        return Err(ResolutionError::UnsortedEnergies.into());
290    }
291    let ext_grid = build_aux_grid(energies, instrument, resonance_data);
292    let (_, layout) = working_grid_layout(energies, ext_grid.as_ref());
293    Ok(layout)
294}
295
296/// Errors from the transmission forward model.
297#[derive(Debug)]
298pub enum TransmissionError {
299    /// The energy grid is not sorted or has a length mismatch with data.
300    Resolution(ResolutionError),
301    /// Doppler broadening parameter validation failed.
302    Doppler(DopplerParamsError),
303    /// Doppler broadening input validation failed (e.g. length mismatch).
304    DopplerBroadening(crate::doppler::DopplerError),
305    /// Computation was cancelled via the cancel token.
306    Cancelled,
307    /// Input array mismatch (e.g. cross-sections vs thicknesses length).
308    InputMismatch(String),
309}
310
311impl fmt::Display for TransmissionError {
312    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
313        match self {
314            Self::Resolution(e) => write!(f, "resolution broadening error: {}", e),
315            Self::Doppler(e) => write!(f, "Doppler parameter error: {}", e),
316            Self::DopplerBroadening(e) => write!(f, "Doppler broadening error: {}", e),
317            Self::Cancelled => write!(f, "computation cancelled"),
318            Self::InputMismatch(msg) => write!(f, "input mismatch: {}", msg),
319        }
320    }
321}
322
323impl std::error::Error for TransmissionError {
324    fn source(&self) -> Option<&(dyn std::error::Error + 'static)> {
325        match self {
326            Self::Resolution(e) => Some(e),
327            Self::Doppler(e) => Some(e),
328            Self::DopplerBroadening(e) => Some(e),
329            Self::Cancelled => None,
330            Self::InputMismatch(_) => None,
331        }
332    }
333}
334
335impl From<ResolutionError> for TransmissionError {
336    fn from(e: ResolutionError) -> Self {
337        Self::Resolution(e)
338    }
339}
340
341impl From<DopplerParamsError> for TransmissionError {
342    fn from(e: DopplerParamsError) -> Self {
343        Self::Doppler(e)
344    }
345}
346
347impl From<crate::doppler::DopplerError> for TransmissionError {
348    fn from(e: crate::doppler::DopplerError) -> Self {
349        Self::DopplerBroadening(e)
350    }
351}
352
353/// Broadened cross-sections and their temperature derivative.
354///
355/// `xs[k][e]` is the Doppler+resolution-broadened cross-section for isotope
356/// `k` at energy index `e`; `dxs_dt[k][e]` is the analytical derivative
357/// with respect to temperature.
358pub type BroadenedXsWithDerivative = (Vec<Vec<f64>>, Vec<Vec<f64>>);
359
360/// The working energy grid used for broadening, plus the map back to the
361/// data grid.
362///
363/// `forward_model` and the Beer-Lambert-aware transmission pipeline run
364/// Doppler/Beer-Lambert/resolution on a *working* grid — the auxiliary
365/// extended grid (boundary extension + resonance fine-structure) when a
366/// Gaussian resolution function is active, otherwise the data grid itself —
367/// and extract the data points LAST.  Fitting models that cache broadened σ
368/// for reuse across LM steps need this layout so they can reproduce the same
369/// "broaden-on-working-grid, extract-last" ordering (issue #608): the LM fit's
370/// cached / precomputed paths previously collapsed σ to the coarse data grid
371/// *before* resolution broadening, degrading the convolution near grid edges
372/// and around narrow resonances relative to [`forward_model`].
373///
374/// **Tabulated resolution has no auxiliary grid.**  `build_aux_grid` only
375/// extends the grid for [`ResolutionFunction::Gaussian`]; for tabulated
376/// kernels (and when no instrument is present) `energies` equals the data
377/// grid and `data_indices` is the identity `0..n`.  Callers can therefore
378/// continue to use a data-grid [`crate::resolution::ResolutionPlan`] for
379/// tabulated kernels — the working grid and the data grid coincide, so the
380/// plan's grid-identity check still passes.
381#[derive(Debug, Clone)]
382pub struct WorkingGridLayout {
383    /// Working-grid energies (eV, ascending).  Equals the input data grid
384    /// when no auxiliary grid was built (tabulated/no resolution).
385    pub energies: Vec<f64>,
386    /// `data_indices[i]` is the index of data energy `i` within
387    /// [`Self::energies`].  Identity (`0..n`) when no auxiliary grid.
388    pub data_indices: Vec<usize>,
389}
390
391impl WorkingGridLayout {
392    /// `true` when the working grid is the data grid itself (no auxiliary
393    /// extension was built).  In that case [`Self::extract`] is a no-op clone.
394    pub fn is_identity(&self) -> bool {
395        self.data_indices.len() == self.energies.len()
396            && self
397                .data_indices
398                .iter()
399                .enumerate()
400                .all(|(i, &idx)| i == idx)
401    }
402
403    /// Extract the data-grid points from a working-grid spectrum.
404    pub fn extract(&self, working: &[f64]) -> Vec<f64> {
405        self.data_indices.iter().map(|&i| working[i]).collect()
406    }
407}
408
409/// Working-grid Doppler-broadened cross-sections plus the data-grid map.
410///
411/// Resolution broadening is **not** applied (issue #442: resolution is applied
412/// to the total transmission after Beer-Lambert).  The caller applies
413/// Beer-Lambert + resolution on [`WorkingGridLayout::energies`] and extracts
414/// the data points last via [`WorkingGridLayout::extract`].
415pub struct WorkingGridXs {
416    /// Doppler-broadened σ per isotope, on the working grid.
417    pub sigma: Vec<Vec<f64>>,
418    /// Working-grid layout (energies + data-grid index map).
419    pub layout: WorkingGridLayout,
420}
421
422/// Working-grid Doppler-broadened cross-sections and their analytical
423/// temperature derivative, plus the data-grid map.
424///
425/// Like [`WorkingGridXs`] but additionally carries `∂σ/∂T` on the working
426/// grid so the analytical Jacobian can form `−T·Σᵢ nᵢ·rᵢ·∂σᵢ/∂T` on the
427/// working grid, resolution-broaden it there, and extract the data points
428/// last (issue #608).
429pub struct WorkingGridXsWithDerivative {
430    /// Doppler-broadened σ per isotope, on the working grid.
431    pub sigma: Vec<Vec<f64>>,
432    /// `∂σ/∂T` per isotope, on the working grid.
433    pub dsigma_dt: Vec<Vec<f64>>,
434    /// Working-grid layout (energies + data-grid index map).
435    pub layout: WorkingGridLayout,
436}
437
438/// Compute transmission from cross-sections via Beer-Lambert law.
439///
440/// T(E) = exp(-thickness × σ(E))
441///
442/// # Arguments
443/// * `cross_sections` — Total cross-sections in barns at each energy point.
444/// * `thickness` — Areal density in atoms/barn (= number_density × path_length).
445///
446/// # Returns
447/// Transmission values (0 to 1) at each energy point.
448pub fn beer_lambert(cross_sections: &[f64], thickness: f64) -> Vec<f64> {
449    cross_sections
450        .iter()
451        .map(|&sigma| (-thickness * sigma).exp())
452        .collect()
453}
454
455/// Compute transmission for multiple isotopes.
456///
457/// T(E) = exp(-Σᵢ thicknessᵢ × σᵢ(E))
458///
459/// # Arguments
460/// * `cross_sections_per_isotope` — Vec of cross-section arrays, one per isotope.
461///   Each inner slice has the same length as the energy grid.
462/// * `thicknesses` — Areal density (atoms/barn) for each isotope.
463///
464/// # Returns
465/// Combined transmission values at each energy point.
466pub fn beer_lambert_multi(
467    cross_sections_per_isotope: &[&[f64]],
468    thicknesses: &[f64],
469) -> Result<Vec<f64>, TransmissionError> {
470    if cross_sections_per_isotope.len() != thicknesses.len() {
471        return Err(TransmissionError::InputMismatch(format!(
472            "cross_sections_per_isotope length ({}) must match thicknesses length ({})",
473            cross_sections_per_isotope.len(),
474            thicknesses.len()
475        )));
476    }
477    if cross_sections_per_isotope.is_empty() {
478        return Err(TransmissionError::InputMismatch(
479            "cross_sections_per_isotope must not be empty".into(),
480        ));
481    }
482
483    let n_energies = cross_sections_per_isotope[0].len();
484    for (k, sigma) in cross_sections_per_isotope.iter().enumerate() {
485        if sigma.len() != n_energies {
486            return Err(TransmissionError::InputMismatch(format!(
487                "cross_sections_per_isotope[{}] length ({}) must match [0] length ({})",
488                k,
489                sigma.len(),
490                n_energies
491            )));
492        }
493    }
494
495    Ok((0..n_energies)
496        .map(|i| {
497            let total_attenuation: f64 = cross_sections_per_isotope
498                .iter()
499                .zip(thicknesses.iter())
500                .map(|(sigma, &thick)| thick * sigma[i])
501                .sum();
502            (-total_attenuation).exp()
503        })
504        .collect())
505}
506
507/// Errors from `SampleParams` construction.
508#[derive(Debug, PartialEq)]
509pub enum SampleParamsError {
510    /// Temperature must be finite.
511    NonFiniteTemperature(f64),
512    /// Temperature must be non-negative.
513    NegativeTemperature(f64),
514}
515
516impl fmt::Display for SampleParamsError {
517    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
518        match self {
519            Self::NonFiniteTemperature(v) => {
520                write!(f, "temperature must be finite, got {v}")
521            }
522            Self::NegativeTemperature(v) => {
523                write!(f, "temperature must be non-negative, got {v}")
524            }
525        }
526    }
527}
528
529impl std::error::Error for SampleParamsError {}
530
531/// Sample description for the forward model.
532#[derive(Debug, Clone)]
533pub struct SampleParams {
534    /// Temperature in Kelvin (for Doppler broadening).
535    temperature_k: f64,
536    /// Isotope compositions: (resonance data, areal density in atoms/barn).
537    isotopes: Vec<(ResonanceData, f64)>,
538}
539
540impl SampleParams {
541    /// Create validated sample parameters.
542    ///
543    /// # Errors
544    /// Returns `SampleParamsError::NonFiniteTemperature` if `temperature_k` is
545    /// NaN or infinity.
546    /// Returns `SampleParamsError::NegativeTemperature` if `temperature_k < 0.0`.
547    pub fn new(
548        temperature_k: f64,
549        isotopes: Vec<(ResonanceData, f64)>,
550    ) -> Result<Self, SampleParamsError> {
551        if !temperature_k.is_finite() {
552            return Err(SampleParamsError::NonFiniteTemperature(temperature_k));
553        }
554        if temperature_k < 0.0 {
555            return Err(SampleParamsError::NegativeTemperature(temperature_k));
556        }
557        Ok(Self {
558            temperature_k,
559            isotopes,
560        })
561    }
562
563    /// Returns the sample temperature in Kelvin.
564    #[must_use]
565    pub fn temperature_k(&self) -> f64 {
566        self.temperature_k
567    }
568
569    /// Returns the isotope compositions: (resonance data, areal density).
570    #[must_use]
571    pub fn isotopes(&self) -> &[(ResonanceData, f64)] {
572        &self.isotopes
573    }
574}
575
576/// Optional instrument resolution parameters.
577#[derive(Debug, Clone)]
578pub struct InstrumentParams {
579    /// Resolution broadening function (Gaussian or tabulated).
580    pub resolution: ResolutionFunction,
581}
582
583/// Compute a complete theoretical transmission spectrum.
584///
585/// This is the main forward model that chains:
586///   ENDF parameters → cross-sections → Doppler broadening → resolution → transmission
587///
588/// # Arguments
589/// * `energies` — Energy grid in eV (sorted ascending).
590/// * `sample` — Sample parameters (isotopes with areal densities, temperature).
591/// * `instrument` — Optional instrument parameters (resolution broadening).
592///
593/// # Returns
594/// Theoretical transmission spectrum on the energy grid.
595///
596/// # Errors
597/// * [`TransmissionError::Resolution`] — if resolution broadening is
598///   enabled (`instrument` is `Some`) and `energies` is not sorted ascending.
599/// * [`TransmissionError::Doppler`] — if Doppler broadening is enabled
600///   (`temperature_k > 0.0`) and `DopplerParams` validation fails
601///   (e.g., non-positive or non-finite AWR).
602///
603/// **Note**: isotopes with thickness <= 0.0 are silently skipped
604/// (they contribute zero attenuation). This allows callers to include
605/// inactive isotopes in `SampleParams` without causing errors.
606pub fn forward_model(
607    energies: &[f64],
608    sample: &SampleParams,
609    instrument: Option<&InstrumentParams>,
610) -> Result<Vec<f64>, TransmissionError> {
611    let n = energies.len();
612    if n == 0 {
613        return Ok(vec![]);
614    }
615
616    // Validate energy grid once before the per-isotope loop so that
617    // resolution broadening can use the presorted (unchecked) path,
618    // avoiding redundant O(N) sort checks per isotope.
619    if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
620        return Err(ResolutionError::UnsortedEnergies.into());
621    }
622
623    // Build auxiliary grid with boundary extension + resonance fine-structure.
624    // Collect references to avoid cloning full ResonanceData structs.
625    // SAMMY Ref: dat/mdat4.f90 Escale, Fspken, Add_Pnts
626    let active_rd: Vec<&ResonanceData> = sample
627        .isotopes()
628        .iter()
629        .filter(|(_, t)| *t > 0.0)
630        .map(|(rd, _)| rd)
631        .collect();
632    let ext_grid = build_aux_grid(energies, instrument, &active_rd);
633
634    // Compute Doppler-broadened cross-sections for all isotopes in parallel.
635    // Resolution is NOT applied here — it must be applied after Beer-Lambert
636    // on the total transmission.
637    //
638    // SAMMY Ref: DopplerAndResolutionBroadener.cpp — resolution broadening is
639    // applied to T(E), not to σ(E).  Due to Jensen's inequality (exp is
640    // convex), broadening σ before the nonlinear Beer-Lambert systematically
641    // overestimates effective cross-sections at resonance peaks.
642    //
643    // Correct pipeline:
644    //   1. Per-isotope: σ → Doppler → σ_D   (on working grid)
645    //   2. Accumulate:  attenuation = Σᵢ nᵢ·σ_{D,i}
646    //   3. Beer-Lambert: T = exp(−attenuation)
647    //   4. Resolution:  T_broad = R ⊗ T     (on working grid)
648    //   5. Extract at data positions
649    //
650    // The working grid is the extended grid (with boundary+fine-structure
651    // points) when available, otherwise the data grid.
652
653    // Determine working grid: extended grid for resolution boundary handling,
654    // or the data grid when no extension was needed.
655    let (work_energies, work_len): (&[f64], usize) = if let Some((ref ext_e, _)) = ext_grid {
656        (ext_e.as_slice(), ext_e.len())
657    } else {
658        (energies, n)
659    };
660
661    let doppler_xs: Result<Vec<(Vec<f64>, f64)>, TransmissionError> = sample
662        .isotopes()
663        .par_iter()
664        .filter(|(_, thickness)| *thickness > 0.0)
665        .map(|(res_data, thickness)| {
666            let unbroadened: Vec<f64> = work_energies
667                .iter()
668                .map(|&e| reich_moore::cross_sections_at_energy(res_data, e).total)
669                .collect();
670            let after_doppler = if sample.temperature_k() > 0.0 {
671                let params = DopplerParams::new(sample.temperature_k(), res_data.awr)?;
672                doppler::doppler_broaden(work_energies, &unbroadened, &params)?
673            } else {
674                unbroadened
675            };
676            Ok((after_doppler, *thickness))
677        })
678        .collect();
679    let doppler_xs = doppler_xs?;
680
681    // 4. Accumulate total attenuation: Σᵢ thicknessᵢ × σ_{D,i}(E)
682    let mut total_attenuation = vec![0.0f64; work_len];
683    for (xs, thickness) in &doppler_xs {
684        for i in 0..work_len {
685            total_attenuation[i] += thickness * xs[i];
686        }
687    }
688
689    // 5. Beer-Lambert: T = exp(−attenuation)
690    let transmission: Vec<f64> = total_attenuation.iter().map(|&att| (-att).exp()).collect();
691
692    // 6. Resolution broadening on total transmission, then extract at data positions.
693    if let Some(inst) = instrument {
694        let t_broadened =
695            resolution::apply_resolution_presorted(work_energies, &transmission, &inst.resolution);
696        if let Some((_, ref data_indices)) = ext_grid {
697            Ok(data_indices.iter().map(|&i| t_broadened[i]).collect())
698        } else {
699            Ok(t_broadened)
700        }
701    } else {
702        Ok(transmission)
703    }
704}
705
706/// Compute Doppler-broadened cross-sections for each isotope.
707///
708/// Returns **Doppler-only** cross-sections.  Resolution broadening is NOT
709/// applied here because it must be applied after Beer-Lambert on the total
710/// transmission for physically correct results (issue #442).
711///
712/// When `instrument` is `Some`, the auxiliary extended grid is still
713/// constructed for Doppler boundary accuracy, but the resolution
714/// convolution is not performed.
715///
716/// This is the expensive physics step; cache the result once before
717/// fitting many pixels on the same isotopes and energy grid.
718///
719/// # Arguments
720/// * `energies`        — Energy grid in eV (sorted ascending).
721/// * `resonance_data`  — Resonance parameters for each isotope.
722/// * `temperature_k`   — Sample temperature for Doppler broadening.
723/// * `instrument`      — Optional instrument resolution parameters.
724///   Used only for auxiliary grid construction (Doppler boundary accuracy).
725///   Resolution broadening is NOT applied.
726/// * `cancel`          — Optional cancellation token.  Cancellation is checked
727///   at the start of each isotope's parallel task; in-flight tasks run to
728///   completion (consistent with the rayon pattern in `spatial.rs`).
729///
730/// # Returns
731/// One cross-section vector per isotope on success.
732///
733/// # Errors
734/// * [`TransmissionError::Cancelled`] — if the `cancel` flag was observed
735///   during parallel execution (either before an isotope started or after
736///   all tasks completed).
737/// * [`TransmissionError::Resolution`] — if `instrument` is `Some` and
738///   `energies` is not sorted ascending.
739/// * [`TransmissionError::Doppler`] — if Doppler broadening is enabled
740///   (`temperature_k > 0.0`) and `DopplerParams` validation fails
741///   (e.g., non-positive or non-finite AWR).
742pub fn broadened_cross_sections(
743    energies: &[f64],
744    resonance_data: &[ResonanceData],
745    temperature_k: f64,
746    instrument: Option<&InstrumentParams>,
747    cancel: Option<&AtomicBool>,
748) -> Result<Vec<Vec<f64>>, TransmissionError> {
749    // Delegate to the working-grid variant and extract the data points.
750    // Doppler runs on the working (aux) grid for boundary accuracy.
751    let WorkingGridXs { sigma, layout } = broadened_cross_sections_on_working_grid(
752        energies,
753        resonance_data,
754        temperature_k,
755        instrument,
756        cancel,
757    )?;
758    Ok(sigma.iter().map(|s| layout.extract(s)).collect())
759}
760
761/// Like [`broadened_cross_sections`] but returns the Doppler-broadened σ on the
762/// **working grid** (auxiliary extended grid when Gaussian resolution is active,
763/// else the data grid) together with the [`WorkingGridLayout`].
764///
765/// The spatial production pipeline (`spatial_map_typed`) uses this to pre-store
766/// σ on the working grid so each per-pixel `PrecomputedTransmissionModel`
767/// applies Beer-Lambert + resolution on the working grid and extracts the data
768/// points last — matching [`forward_model`] (issue #608).  Resolution is NOT
769/// applied (issue #442).
770pub fn broadened_cross_sections_on_working_grid(
771    energies: &[f64],
772    resonance_data: &[ResonanceData],
773    temperature_k: f64,
774    instrument: Option<&InstrumentParams>,
775    cancel: Option<&AtomicBool>,
776) -> Result<WorkingGridXs, TransmissionError> {
777    // Validate energy grid once before the per-isotope loop.
778    if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
779        return Err(ResolutionError::UnsortedEnergies.into());
780    }
781
782    // Build auxiliary grid with boundary extension + resonance fine-structure.
783    // SAMMY extends the energy grid beyond the data range and adds dense points
784    // around narrow resonances so the broadening convolution integrals have
785    // adequate quadrature points.
786    // SAMMY Ref: dat/mdat4.f90 Escale+Fspken+Add_Pnts, dat/mdata.f90 Vqcon
787    let rd_refs: Vec<&ResonanceData> = resonance_data.iter().collect();
788    let ext_grid = build_aux_grid(energies, instrument, &rd_refs);
789    let (work_energies, layout) = working_grid_layout(energies, ext_grid.as_ref());
790
791    // Parallelize across isotopes — Doppler broadening for each isotope is
792    // independent.  Resolution is NOT applied here (issue #442: resolution
793    // must be applied after Beer-Lambert on total transmission).  σ is returned
794    // on the working grid WITHOUT extracting the data points.
795    // Cancellation is checked per-isotope inside the parallel map.
796    let result: Result<Vec<Vec<f64>>, TransmissionError> = resonance_data
797        .par_iter()
798        .map(|rd| {
799            // Check cancellation before starting this isotope.
800            if cancel.is_some_and(|c| c.load(Ordering::Relaxed)) {
801                return Err(TransmissionError::Cancelled);
802            }
803
804            let unbroadened: Vec<f64> = work_energies
805                .iter()
806                .map(|&e| reich_moore::cross_sections_at_energy(rd, e).total)
807                .collect();
808            if temperature_k > 0.0 {
809                let params = DopplerParams::new(temperature_k, rd.awr)?;
810                doppler::doppler_broaden(work_energies, &unbroadened, &params).map_err(Into::into)
811            } else {
812                Ok(unbroadened)
813            }
814        })
815        .collect();
816
817    // Final cancellation check: if cancel was set during parallel execution,
818    // some tasks may have completed before observing it.
819    if cancel.is_some_and(|c| c.load(Ordering::Relaxed)) {
820        return Err(TransmissionError::Cancelled);
821    }
822
823    Ok(WorkingGridXs {
824        sigma: result?,
825        layout,
826    })
827}
828
829/// Compute Doppler+resolution-broadened cross-sections using SAMMY's
830/// Beer-Lambert-aware pipeline for transmission data.
831///
832/// For transmission data, SAMMY applies resolution broadening to the
833/// transmission T = exp(-nd×σ_D) rather than to σ_D directly.  Due to
834/// Jensen's inequality (the exponential is convex), direct σ broadening
835/// overestimates the effective cross section at resonance peaks.  This
836/// function implements SAMMY's correct pipeline:
837///
838/// 1. Evaluate unbroadened σ on extended grid
839/// 2. Doppler-broaden σ → σ_D
840/// 3. Convert to transmission: T = exp(-nd × σ_D)
841/// 4. Resolution-broaden T → T_broadened
842/// 5. Convert back: σ_eff = -ln(T_broadened) / nd
843///
844/// SAMMY Ref: DopplerAndResolutionBroadener.cpp — resolution broadening is
845/// applied after Beer-Lambert conversion in the SAMMY pipeline.
846///
847/// # Arguments
848/// * `energies`             — Energy grid in eV (sorted ascending).
849/// * `resonance_data`       — Resonance parameters for each isotope.
850/// * `temperature_k`        — Sample temperature for Doppler broadening.
851/// * `instrument`           — Instrument resolution parameters.
852/// * `thickness_atoms_barn`  — Sample thickness n×d (atoms/barn).  Must be > 0.
853/// * `cancel`               — Optional cancellation token.
854///
855/// # Returns
856/// One effective cross-section vector per isotope on success.
857pub fn broadened_cross_sections_for_transmission(
858    energies: &[f64],
859    resonance_data: &[ResonanceData],
860    temperature_k: f64,
861    instrument: &InstrumentParams,
862    thickness_atoms_barn: f64,
863    cancel: Option<&AtomicBool>,
864) -> Result<Vec<Vec<f64>>, TransmissionError> {
865    // The Beer-Lambert→σ_eff inversion (step 5: σ_eff = −ln(T)/nd) divides by
866    // the thickness, so a non-positive or non-finite nd would yield ±∞/NaN
867    // cross-sections.  The docstring contract requires nd > 0; enforce it
868    // up-front before any work, matching the module's `InputMismatch` style.
869    if !thickness_atoms_barn.is_finite() || thickness_atoms_barn <= 0.0 {
870        return Err(TransmissionError::InputMismatch(format!(
871            "thickness_atoms_barn must be finite and > 0, got {thickness_atoms_barn}"
872        )));
873    }
874    if !energies.windows(2).all(|w| w[0] <= w[1]) {
875        return Err(ResolutionError::UnsortedEnergies.into());
876    }
877
878    let rd_refs: Vec<&ResonanceData> = resonance_data.iter().collect();
879    let ext_grid = build_aux_grid(energies, Some(instrument), &rd_refs);
880    let nd = thickness_atoms_barn;
881
882    let result: Result<Vec<Vec<f64>>, TransmissionError> = resonance_data
883        .par_iter()
884        .map(|rd| {
885            if cancel.is_some_and(|c| c.load(Ordering::Relaxed)) {
886                return Err(TransmissionError::Cancelled);
887            }
888
889            let sigma_eff = if let Some((ref ext_energies, ref data_indices)) = ext_grid {
890                // Extended grid available: evaluate the full pipeline on the
891                // extended grid and extract at data positions.
892
893                // 1. Unbroadened cross sections on extended grid.
894                let unbroadened: Vec<f64> = ext_energies
895                    .iter()
896                    .map(|&e| reich_moore::cross_sections_at_energy(rd, e).total)
897                    .collect();
898
899                // 2. Doppler broadening.
900                let after_doppler = if temperature_k > 0.0 {
901                    let params = DopplerParams::new(temperature_k, rd.awr)?;
902                    doppler::doppler_broaden(ext_energies, &unbroadened, &params)?
903                } else {
904                    unbroadened
905                };
906
907                // 3. Convert to transmission: T = exp(-nd × σ_D).
908                let transmission: Vec<f64> = after_doppler
909                    .iter()
910                    .map(|&sigma| (-nd * sigma).exp())
911                    .collect();
912
913                // 4. Resolution-broaden T.
914                let t_broadened = resolution::apply_resolution_presorted(
915                    ext_energies,
916                    &transmission,
917                    &instrument.resolution,
918                );
919
920                // 5. Convert back to effective σ: σ_eff = -ln(T_broad) / nd.
921                data_indices
922                    .iter()
923                    .map(|&i| {
924                        let t = t_broadened[i].clamp(1e-30, 1.0);
925                        -t.ln() / nd
926                    })
927                    .collect()
928            } else {
929                // No extended grid (e.g. tabulated resolution with no aux grid):
930                // Doppler on data grid, Beer-Lambert, resolution on data grid.
931                let unbroadened: Vec<f64> = energies
932                    .iter()
933                    .map(|&e| reich_moore::cross_sections_at_energy(rd, e).total)
934                    .collect();
935
936                let after_doppler = if temperature_k > 0.0 {
937                    let params = DopplerParams::new(temperature_k, rd.awr)?;
938                    doppler::doppler_broaden(energies, &unbroadened, &params)?
939                } else {
940                    unbroadened
941                };
942
943                let transmission: Vec<f64> = after_doppler
944                    .iter()
945                    .map(|&sigma| (-nd * sigma).exp())
946                    .collect();
947
948                let t_broadened = resolution::apply_resolution_presorted(
949                    energies,
950                    &transmission,
951                    &instrument.resolution,
952                );
953
954                t_broadened
955                    .iter()
956                    .map(|&t| {
957                        let t_clamped = t.clamp(1e-30, 1.0);
958                        -t_clamped.ln() / nd
959                    })
960                    .collect()
961            };
962
963            Ok(sigma_eff)
964        })
965        .collect();
966
967    if cancel.is_some_and(|c| c.load(Ordering::Relaxed)) {
968        return Err(TransmissionError::Cancelled);
969    }
970
971    result
972}
973
974/// Compute unbroadened (raw Reich-Moore) cross-sections for each isotope.
975///
976/// This is the temperature-independent first step of the forward model.
977/// The result can be cached and reused across multiple temperature evaluations
978/// (e.g., during LM iterations where temperature is a free parameter).
979///
980/// # Returns
981/// One total cross-section vector per isotope: `result[k][e]` is the
982/// unbroadened total cross-section (barns) for isotope `k` at energy `e`.
983pub fn unbroadened_cross_sections(
984    energies: &[f64],
985    resonance_data: &[ResonanceData],
986    cancel: Option<&AtomicBool>,
987) -> Result<Vec<Vec<f64>>, TransmissionError> {
988    let result: Result<Vec<Vec<f64>>, TransmissionError> = resonance_data
989        .par_iter()
990        .map(|rd| {
991            if cancel.is_some_and(|c| c.load(Ordering::Relaxed)) {
992                return Err(TransmissionError::Cancelled);
993            }
994            let xs: Vec<f64> = energies
995                .iter()
996                .map(|&e| reich_moore::cross_sections_at_energy(rd, e).total)
997                .collect();
998            Ok(xs)
999        })
1000        .collect();
1001
1002    if cancel.is_some_and(|c| c.load(Ordering::Relaxed)) {
1003        return Err(TransmissionError::Cancelled);
1004    }
1005    result
1006}
1007
1008/// Compute Doppler-broadened cross-sections from precomputed unbroadened
1009/// cross-sections.
1010///
1011/// Returns **Doppler-only** cross-sections.  Resolution broadening is NOT
1012/// applied (issue #442: must be applied after Beer-Lambert on total T).
1013///
1014/// Like [`broadened_cross_sections`] but skips the expensive Reich-Moore
1015/// calculation (step 1). Use [`unbroadened_cross_sections`] to compute
1016/// `base_xs` once, then call this function repeatedly with different
1017/// temperatures.
1018pub fn broadened_cross_sections_from_base(
1019    energies: &[f64],
1020    base_xs: &[Vec<f64>],
1021    resonance_data: &[ResonanceData],
1022    temperature_k: f64,
1023    instrument: Option<&InstrumentParams>,
1024) -> Result<Vec<Vec<f64>>, TransmissionError> {
1025    // Delegate to the working-grid variant and extract the data points.
1026    // Doppler runs on the working (aux) grid for boundary accuracy; resolution
1027    // is NOT applied here (issue #442).  Extracting the data points here keeps
1028    // this function's contract (data-grid σ) unchanged.
1029    let WorkingGridXs { sigma, layout } = broadened_cross_sections_from_base_on_working_grid(
1030        energies,
1031        base_xs,
1032        resonance_data,
1033        temperature_k,
1034        instrument,
1035    )?;
1036    Ok(sigma.iter().map(|s| layout.extract(s)).collect())
1037}
1038
1039/// Like [`broadened_cross_sections_from_base`] but returns the Doppler-broadened
1040/// σ on the **working grid** (auxiliary extended grid when Gaussian resolution
1041/// is active, else the data grid) together with the [`WorkingGridLayout`].
1042///
1043/// This is the building block the LM fit's cached / precomputed paths use to
1044/// apply resolution broadening on the working grid and extract the data points
1045/// last — matching [`forward_model`] (issue #608).  Resolution is NOT applied
1046/// (issue #442).
1047pub fn broadened_cross_sections_from_base_on_working_grid(
1048    energies: &[f64],
1049    base_xs: &[Vec<f64>],
1050    resonance_data: &[ResonanceData],
1051    temperature_k: f64,
1052    instrument: Option<&InstrumentParams>,
1053) -> Result<WorkingGridXs, TransmissionError> {
1054    validate_base_xs(energies, base_xs, resonance_data)?;
1055    if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
1056        return Err(ResolutionError::UnsortedEnergies.into());
1057    }
1058
1059    // Build auxiliary grid with boundary extension + resonance fine-structure.
1060    // base_xs is on the data grid; we extend it to the aux grid by evaluating
1061    // cross-sections at the auxiliary-only points (cheap: only the few hundred
1062    // extra points, not the full grid).
1063    // SAMMY Ref: dat/mdat4.f90 Escale+Fspken+Add_Pnts
1064    let rd_refs: Vec<&ResonanceData> = resonance_data.iter().collect();
1065    let ext_grid = build_aux_grid(energies, instrument, &rd_refs);
1066    let is_data_point = data_point_mask(ext_grid.as_ref());
1067    let (work_energies, layout) = working_grid_layout(energies, ext_grid.as_ref());
1068
1069    // Resolution is NOT applied (issue #442).  Doppler broadening on the
1070    // working grid, returned WITHOUT extracting the data points.
1071    let sigma: Result<Vec<Vec<f64>>, TransmissionError> = base_xs
1072        .par_iter()
1073        .zip(resonance_data.par_iter())
1074        .map(|(xs_raw, rd)| {
1075            let xs_work = if let Some((ref ext_energies, ref data_indices)) = ext_grid {
1076                let mask = is_data_point.as_ref().unwrap();
1077                build_extended_xs_from_base(ext_energies, data_indices, mask, xs_raw, rd)
1078            } else {
1079                xs_raw.clone()
1080            };
1081            if temperature_k > 0.0 {
1082                let params = DopplerParams::new(temperature_k, rd.awr)?;
1083                doppler::doppler_broaden(work_energies, &xs_work, &params).map_err(Into::into)
1084            } else {
1085                Ok(xs_work)
1086            }
1087        })
1088        .collect();
1089
1090    Ok(WorkingGridXs {
1091        sigma: sigma?,
1092        layout,
1093    })
1094}
1095
1096/// Compute Doppler-broadened cross-sections and their **analytical**
1097/// temperature derivative from precomputed unbroadened cross-sections.
1098///
1099/// Returns **Doppler-only** cross-sections and derivatives.  Resolution
1100/// broadening is NOT applied (issue #442: must be applied after
1101/// Beer-Lambert on total T).
1102///
1103/// Uses `doppler_broaden_with_derivative` for exact ∂σ/∂T in a single pass
1104/// (1× broadening), replacing the 3× FD approach.
1105///
1106/// Returns `BroadenedXsWithDerivative`: `(sigma_k, dsigma_k_dT)`.
1107pub fn broadened_cross_sections_with_analytical_derivative_from_base(
1108    energies: &[f64],
1109    base_xs: &[Vec<f64>],
1110    resonance_data: &[ResonanceData],
1111    temperature_k: f64,
1112    instrument: Option<&InstrumentParams>,
1113) -> Result<BroadenedXsWithDerivative, TransmissionError> {
1114    // Delegate to the working-grid variant and extract the data points for
1115    // both σ and ∂σ/∂T.  Keeps this function's data-grid contract unchanged.
1116    let WorkingGridXsWithDerivative {
1117        sigma,
1118        dsigma_dt,
1119        layout,
1120    } = broadened_cross_sections_with_analytical_derivative_from_base_on_working_grid(
1121        energies,
1122        base_xs,
1123        resonance_data,
1124        temperature_k,
1125        instrument,
1126    )?;
1127    let xs_all = sigma.iter().map(|s| layout.extract(s)).collect();
1128    let dxs_all = dsigma_dt.iter().map(|d| layout.extract(d)).collect();
1129    Ok((xs_all, dxs_all))
1130}
1131
1132/// Like [`broadened_cross_sections_with_analytical_derivative_from_base`] but
1133/// returns σ and ∂σ/∂T on the **working grid** (auxiliary extended grid when
1134/// Gaussian resolution is active, else the data grid) together with the
1135/// [`WorkingGridLayout`].
1136///
1137/// The analytical Jacobian uses this so it can form `−T·Σᵢ nᵢ·rᵢ·∂σᵢ/∂T` on
1138/// the working grid, resolution-broaden it there, and extract the data points
1139/// last (issue #608).  Resolution is NOT applied (issue #442).
1140pub fn broadened_cross_sections_with_analytical_derivative_from_base_on_working_grid(
1141    energies: &[f64],
1142    base_xs: &[Vec<f64>],
1143    resonance_data: &[ResonanceData],
1144    temperature_k: f64,
1145    instrument: Option<&InstrumentParams>,
1146) -> Result<WorkingGridXsWithDerivative, TransmissionError> {
1147    validate_base_xs(energies, base_xs, resonance_data)?;
1148    if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
1149        return Err(ResolutionError::UnsortedEnergies.into());
1150    }
1151
1152    // Build auxiliary grid (same as broadened_cross_sections_from_base).
1153    let rd_refs: Vec<&ResonanceData> = resonance_data.iter().collect();
1154    let ext_grid = build_aux_grid(energies, instrument, &rd_refs);
1155    let is_data_point = data_point_mask(ext_grid.as_ref());
1156    let (work_energies, layout) = working_grid_layout(energies, ext_grid.as_ref());
1157
1158    // Per-isotope: Doppler broaden with analytical derivative on the WORKING
1159    // grid, returned WITHOUT extracting the data points.
1160    // Resolution is NOT applied (issue #442).
1161    type IsotopeXsDxs = Result<(Vec<f64>, Vec<f64>), TransmissionError>;
1162    let results: Vec<IsotopeXsDxs> = base_xs
1163        .par_iter()
1164        .zip(resonance_data.par_iter())
1165        .map(|(xs_raw, rd)| {
1166            let xs_work = if let Some((ref ext_energies, ref data_indices)) = ext_grid {
1167                let mask = is_data_point.as_ref().unwrap();
1168                build_extended_xs_from_base(ext_energies, data_indices, mask, xs_raw, rd)
1169            } else {
1170                xs_raw.clone()
1171            };
1172            if temperature_k > 0.0 {
1173                let params = DopplerParams::new(temperature_k, rd.awr)?;
1174                doppler::doppler_broaden_with_derivative(work_energies, &xs_work, &params)
1175                    .map_err(Into::into)
1176            } else {
1177                let zeros = vec![0.0; work_energies.len()];
1178                Ok((xs_work, zeros))
1179            }
1180        })
1181        .collect();
1182
1183    // Separate into (sigma, dsigma_dt) on the working grid.
1184    let mut sigma = Vec::with_capacity(base_xs.len());
1185    let mut dsigma_dt = Vec::with_capacity(base_xs.len());
1186    for r in results {
1187        let (xs, dxs) = r?;
1188        sigma.push(xs);
1189        dsigma_dt.push(dxs);
1190    }
1191    Ok(WorkingGridXsWithDerivative {
1192        sigma,
1193        dsigma_dt,
1194        layout,
1195    })
1196}
1197
1198/// Compute a transmission spectrum from precomputed unbroadened cross-sections.
1199///
1200/// Applies Doppler broadening and Beer-Lambert using cached base XS, then
1201/// resolution broadening on the total transmission (issue #442).  This skips
1202/// the expensive Reich-Moore calculation, making it suitable for use inside
1203/// `TransmissionFitModel::evaluate()` when temperature is a free parameter.
1204///
1205/// The full pipeline runs on the **auxiliary extended grid** (boundary
1206/// extension + resonance fine-structure points) so that Doppler AND resolution
1207/// broadening have adequate quadrature support and correct boundary handling;
1208/// the data-grid points are extracted LAST.  This mirrors the non-cached
1209/// [`forward_model`] and [`broadened_cross_sections_for_transmission`] exactly
1210/// — previously this cached path collapsed σ to the coarse data grid before
1211/// resolution broadening, degrading the convolution near the grid edges and
1212/// around narrow resonances.
1213///
1214/// Pipeline:
1215///   1. Per isotope: build σ on the extended grid (cached data-grid values +
1216///      freshly-evaluated auxiliary points), Doppler-broaden on the extended grid
1217///   2. Accumulate total attenuation on the extended grid: Σᵢ nᵢ·σ_{D,i}
1218///   3. Beer-Lambert: T = exp(−attenuation) on the extended grid
1219///   4. Resolution: T_broad = R ⊗ T on the extended grid (when instrument present)
1220///   5. Extract the data-grid points
1221pub fn forward_model_from_base_xs(
1222    energies: &[f64],
1223    base_xs: &[Vec<f64>],
1224    resonance_data: &[ResonanceData],
1225    thicknesses: &[f64],
1226    temperature_k: f64,
1227    instrument: Option<&InstrumentParams>,
1228) -> Result<Vec<f64>, TransmissionError> {
1229    if base_xs.len() != resonance_data.len() || thicknesses.len() != resonance_data.len() {
1230        return Err(TransmissionError::InputMismatch(format!(
1231            "forward_model_from_base_xs: base_xs({})/thicknesses({})/resonance_data({}) length mismatch",
1232            base_xs.len(),
1233            thicknesses.len(),
1234            resonance_data.len(),
1235        )));
1236    }
1237    for (i, row) in base_xs.iter().enumerate() {
1238        if row.len() != energies.len() {
1239            return Err(TransmissionError::InputMismatch(format!(
1240                "base_xs[{i}] has {} energies but expected {}",
1241                row.len(),
1242                energies.len(),
1243            )));
1244        }
1245    }
1246    let n = energies.len();
1247    if n == 0 {
1248        return Ok(vec![]);
1249    }
1250
1251    // Validate the data grid once before the per-isotope loop so resolution can
1252    // use the presorted (unchecked) path, matching forward_model.
1253    if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
1254        return Err(ResolutionError::UnsortedEnergies.into());
1255    }
1256
1257    // Build the auxiliary extended grid (boundary + resonance fine-structure).
1258    // SAMMY Ref: dat/mdat4.f90 Escale+Fspken+Add_Pnts.
1259    let rd_refs: Vec<&ResonanceData> = resonance_data.iter().collect();
1260    let ext_grid = build_aux_grid(energies, instrument, &rd_refs);
1261    let is_data_point: Option<Vec<bool>> = ext_grid.as_ref().map(|(ext_e, di)| {
1262        let mut mask = vec![false; ext_e.len()];
1263        for &idx in di {
1264            mask[idx] = true;
1265        }
1266        mask
1267    });
1268
1269    // Working grid: the extended grid when available, else the data grid.
1270    let (work_energies, work_len): (&[f64], usize) = if let Some((ref ext_e, _)) = ext_grid {
1271        (ext_e.as_slice(), ext_e.len())
1272    } else {
1273        (energies, n)
1274    };
1275
1276    // Step 1: per-isotope Doppler-broadened σ on the WORKING (extended) grid.
1277    // Resolution is NOT applied here (issue #442) — it goes on total T below.
1278    let doppler_xs: Result<Vec<Vec<f64>>, TransmissionError> = base_xs
1279        .par_iter()
1280        .zip(resonance_data.par_iter())
1281        .map(|(xs_raw, rd)| {
1282            let xs_ext = if let Some((ref ext_energies, ref data_indices)) = ext_grid {
1283                let mask = is_data_point.as_ref().unwrap();
1284                build_extended_xs_from_base(ext_energies, data_indices, mask, xs_raw, rd)
1285            } else {
1286                xs_raw.clone()
1287            };
1288            if temperature_k > 0.0 {
1289                let params = DopplerParams::new(temperature_k, rd.awr)?;
1290                doppler::doppler_broaden(work_energies, &xs_ext, &params).map_err(Into::into)
1291            } else {
1292                Ok(xs_ext)
1293            }
1294        })
1295        .collect();
1296    let doppler_xs = doppler_xs?;
1297
1298    // Step 2-3: accumulate attenuation on the working grid, then Beer-Lambert.
1299    let mut total_attenuation = vec![0.0f64; work_len];
1300    for (xs, &thickness) in doppler_xs.iter().zip(thicknesses.iter()) {
1301        if thickness <= 0.0 {
1302            continue;
1303        }
1304        for i in 0..work_len {
1305            total_attenuation[i] += thickness * xs[i];
1306        }
1307    }
1308    let transmission: Vec<f64> = total_attenuation.iter().map(|&att| (-att).exp()).collect();
1309
1310    // Step 4-5: resolution broadening on total transmission (on the working
1311    // grid), then extract the data-grid points last.  When `instrument` is
1312    // `None`, `build_aux_grid` returns `None`, so the working grid IS the data
1313    // grid and `transmission` is returned directly (matching `forward_model`).
1314    if let Some(inst) = instrument {
1315        let t_broadened =
1316            resolution::apply_resolution_presorted(work_energies, &transmission, &inst.resolution);
1317        if let Some((_, ref data_indices)) = ext_grid {
1318            Ok(data_indices.iter().map(|&i| t_broadened[i]).collect())
1319        } else {
1320            Ok(t_broadened)
1321        }
1322    } else {
1323        Ok(transmission)
1324    }
1325}
1326
1327#[cfg(test)]
1328mod tests {
1329    use super::*;
1330    use nereids_core::types::Isotope;
1331    use nereids_endf::resonance::test_support::u238_single_resonance;
1332    use nereids_endf::resonance::{LGroup, Resonance, ResonanceFormalism, ResonanceRange};
1333
1334    /// Issue #608 R2: the spatial pipeline builds the working-grid σ
1335    /// (`broadened_cross_sections_on_working_grid`) and the shared layout
1336    /// (`resolution_working_grid`) via SEPARATE calls — they must produce
1337    /// bit-identical layouts, or per-pixel σ and the shared layout would
1338    /// disagree.  Both route through `build_aux_grid` with the same arguments;
1339    /// this pins that determinism against a future divergence.
1340    #[test]
1341    fn working_grid_layout_matches_across_separate_calls() {
1342        let data = u238_single_resonance();
1343        let energies: Vec<f64> = (0..401).map(|i| 4.0 + (i as f64) * 0.015).collect();
1344        let inst = InstrumentParams {
1345            resolution: crate::resolution::ResolutionFunction::Gaussian(
1346                crate::resolution::ResolutionParams::new(25.0, 0.5, 0.005, 0.0).unwrap(),
1347            ),
1348        };
1349        let layout_a = resolution_working_grid(&energies, Some(&inst), &[&data]).unwrap();
1350        let working = broadened_cross_sections_on_working_grid(
1351            &energies,
1352            std::slice::from_ref(&data),
1353            300.0,
1354            Some(&inst),
1355            None,
1356        )
1357        .unwrap();
1358        let layout_b = working.layout;
1359        assert!(
1360            !layout_a.is_identity(),
1361            "Gaussian resolution should build a non-identity auxiliary grid"
1362        );
1363        assert_eq!(
1364            layout_a.data_indices, layout_b.data_indices,
1365            "data-index maps must match across the two builders"
1366        );
1367        assert_eq!(
1368            layout_a.energies.len(),
1369            layout_b.energies.len(),
1370            "working-grid length must match"
1371        );
1372        for (a, b) in layout_a.energies.iter().zip(layout_b.energies.iter()) {
1373            assert_eq!(
1374                a.to_bits(),
1375                b.to_bits(),
1376                "working-grid energies must be bit-identical"
1377            );
1378        }
1379    }
1380
1381    #[test]
1382    fn test_beer_lambert_zero_thickness() {
1383        let xs = vec![100.0, 200.0, 300.0];
1384        let t = beer_lambert(&xs, 0.0);
1385        assert_eq!(t, vec![1.0, 1.0, 1.0]);
1386    }
1387
1388    #[test]
1389    fn test_beer_lambert_basic() {
1390        // σ = 100 barns, thickness = 0.01 atoms/barn
1391        // T = exp(-1.0) ≈ 0.3679
1392        let xs = vec![100.0];
1393        let t = beer_lambert(&xs, 0.01);
1394        assert!(
1395            (t[0] - (-1.0_f64).exp()).abs() < 1e-10,
1396            "T = {}, expected {}",
1397            t[0],
1398            (-1.0_f64).exp()
1399        );
1400    }
1401
1402    #[test]
1403    fn test_beer_lambert_opaque() {
1404        // Very thick sample: T should be 0 (exp(-1000) underflows)
1405        let xs = vec![1000.0];
1406        let t = beer_lambert(&xs, 1.0);
1407        assert_eq!(t[0], 0.0, "T = {}, expected 0.0", t[0]);
1408    }
1409
1410    #[test]
1411    fn test_beer_lambert_multi_additive() {
1412        // Two isotopes should combine additively in the exponent.
1413        // σ₁ = 100 barns, t₁ = 0.01 → att₁ = 1.0
1414        // σ₂ = 200 barns, t₂ = 0.005 → att₂ = 1.0
1415        // T = exp(-(1.0 + 1.0)) = exp(-2.0)
1416        let xs1 = vec![100.0];
1417        let xs2 = vec![200.0];
1418        let t = beer_lambert_multi(&[&xs1, &xs2], &[0.01, 0.005]).unwrap();
1419        assert!(
1420            (t[0] - (-2.0_f64).exp()).abs() < 1e-10,
1421            "T = {}, expected {}",
1422            t[0],
1423            (-2.0_f64).exp()
1424        );
1425    }
1426
1427    #[test]
1428    fn test_transmission_dip_at_resonance() {
1429        // U-238 has a huge capture resonance at 6.674 eV.
1430        // A thin sample should show a transmission dip there.
1431        let data = u238_single_resonance();
1432        let thickness = 0.001; // atoms/barn (thin)
1433
1434        // Evaluate at a few energies
1435        let energies = [1.0, 3.0, 6.674, 10.0, 20.0];
1436        let xs: Vec<f64> = energies
1437            .iter()
1438            .map(|&e| reich_moore::cross_sections_at_energy(&data, e).total)
1439            .collect();
1440        let trans = beer_lambert(&xs, thickness);
1441
1442        // At 6.674 eV (on resonance), transmission should be much lower
1443        let t_on_res = trans[2];
1444        let t_off_res = trans[0]; // 1 eV, off resonance
1445
1446        assert!(
1447            t_on_res < t_off_res,
1448            "On-resonance T ({}) should be < off-resonance T ({})",
1449            t_on_res,
1450            t_off_res
1451        );
1452
1453        // On-resonance with huge σ (~25000 barns), T ≈ exp(-25) ≈ 0
1454        assert!(
1455            t_on_res < 0.01,
1456            "On-resonance T ({}) should be very small",
1457            t_on_res
1458        );
1459    }
1460
1461    #[test]
1462    fn test_forward_model_no_broadening() {
1463        // Forward model at T=0 with no resolution should give
1464        // the same result as direct Beer-Lambert on unbroadened σ.
1465        let data = u238_single_resonance();
1466        let thickness = 0.001;
1467
1468        let energies: Vec<f64> = (0..101).map(|i| 1.0 + (i as f64) * 0.1).collect();
1469
1470        // Direct calculation
1471        let xs: Vec<f64> = energies
1472            .iter()
1473            .map(|&e| reich_moore::cross_sections_at_energy(&data, e).total)
1474            .collect();
1475        let t_direct = beer_lambert(&xs, thickness);
1476
1477        // Forward model
1478        let sample = SampleParams::new(0.0, vec![(data, thickness)]).unwrap();
1479        let t_forward = forward_model(&energies, &sample, None).unwrap();
1480
1481        for i in 0..energies.len() {
1482            assert!(
1483                (t_direct[i] - t_forward[i]).abs() < 1e-10,
1484                "Mismatch at E={}: direct={}, forward={}",
1485                energies[i],
1486                t_direct[i],
1487                t_forward[i]
1488            );
1489        }
1490    }
1491
1492    #[test]
1493    fn test_forward_model_with_broadening() {
1494        // Forward model with Doppler broadening should smooth out the
1495        // transmission dip, making it wider and shallower.
1496        let data = u238_single_resonance();
1497        let thickness = 0.0001; // Very thin (to avoid total absorption)
1498
1499        let energies: Vec<f64> = (0..401).map(|i| 5.0 + (i as f64) * 0.01).collect();
1500
1501        // Cold (no broadening)
1502        let sample_cold = SampleParams::new(0.0, vec![(data.clone(), thickness)]).unwrap();
1503        let t_cold = forward_model(&energies, &sample_cold, None).unwrap();
1504
1505        // Hot (300 K Doppler)
1506        let sample_hot = SampleParams::new(300.0, vec![(data, thickness)]).unwrap();
1507        let t_hot = forward_model(&energies, &sample_hot, None).unwrap();
1508
1509        // Find minima
1510        let min_cold = t_cold.iter().cloned().fold(f64::MAX, f64::min);
1511        let min_hot = t_hot.iter().cloned().fold(f64::MAX, f64::min);
1512
1513        // Broadened dip should be shallower (higher minimum transmission)
1514        assert!(
1515            min_hot > min_cold,
1516            "Broadened min T ({}) should be > unbroadened min T ({})",
1517            min_hot,
1518            min_cold
1519        );
1520    }
1521
1522    #[test]
1523    fn test_forward_model_multi_isotope() {
1524        // Two isotopes with different resonances should create two dips.
1525        let u238 = u238_single_resonance();
1526
1527        // Create a fictitious second isotope with a resonance at 20 eV
1528        let other = ResonanceData {
1529            isotope: Isotope::new(1, 10).unwrap(),
1530            za: 1010,
1531            awr: 10.0,
1532            ranges: vec![ResonanceRange {
1533                energy_low: 0.0,
1534                energy_high: 100.0,
1535                resolved: true,
1536                formalism: ResonanceFormalism::ReichMoore,
1537                target_spin: 0.0,
1538                scattering_radius: 5.0,
1539                naps: 1,
1540                l_groups: vec![LGroup {
1541                    l: 0,
1542                    awr: 10.0,
1543                    apl: 5.0,
1544                    qx: 0.0,
1545                    lrx: 0,
1546                    resonances: vec![Resonance {
1547                        energy: 20.0,
1548                        j: 0.5,
1549                        gn: 0.1,
1550                        gg: 0.05,
1551                        gfa: 0.0,
1552                        gfb: 0.0,
1553                    }],
1554                }],
1555                rml: None,
1556                urr: None,
1557                ap_table: None,
1558                r_external: vec![],
1559            }],
1560        };
1561
1562        let energies: Vec<f64> = (0..301).map(|i| 1.0 + (i as f64) * 0.1).collect();
1563
1564        let sample = SampleParams::new(0.0, vec![(u238, 0.0001), (other, 0.0001)]).unwrap();
1565        let t = forward_model(&energies, &sample, None).unwrap();
1566
1567        // Find the transmission near 6.674 eV (U-238 resonance)
1568        let idx_u238 = energies
1569            .iter()
1570            .position(|&e| (e - 6.7).abs() < 0.05)
1571            .unwrap();
1572        // Find the transmission near 20 eV (other resonance)
1573        let idx_other = energies
1574            .iter()
1575            .position(|&e| (e - 20.0).abs() < 0.05)
1576            .unwrap();
1577        // Off-resonance
1578        let idx_off = energies
1579            .iter()
1580            .position(|&e| (e - 15.0).abs() < 0.05)
1581            .unwrap();
1582
1583        // Both dips should be visible
1584        assert!(
1585            t[idx_u238] < t[idx_off],
1586            "U-238 dip at 6.7 eV: T={}, off-res: T={}",
1587            t[idx_u238],
1588            t[idx_off]
1589        );
1590        assert!(
1591            t[idx_other] < t[idx_off],
1592            "Other dip at 20 eV: T={}, off-res: T={}",
1593            t[idx_other],
1594            t[idx_off]
1595        );
1596    }
1597
1598    #[test]
1599    fn test_broadened_xs_analytical_derivative() {
1600        // Verify analytical ∂σ/∂T against a manual FD at a larger step (1 K)
1601        // and check they agree to reasonable tolerance.
1602        let data = u238_single_resonance();
1603        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
1604        let temperature = 300.0;
1605
1606        let base_xs =
1607            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
1608        let (xs, dxs_dt) = broadened_cross_sections_with_analytical_derivative_from_base(
1609            &energies,
1610            &base_xs,
1611            std::slice::from_ref(&data),
1612            temperature,
1613            None,
1614        )
1615        .unwrap();
1616
1617        // Basic shape checks
1618        assert_eq!(xs.len(), 1, "one isotope");
1619        assert_eq!(dxs_dt.len(), 1, "one isotope derivative");
1620        assert_eq!(xs[0].len(), energies.len());
1621        assert_eq!(dxs_dt[0].len(), energies.len());
1622
1623        // The derivative should be non-zero near the resonance at 6.674 eV
1624        // where Doppler broadening has a strong effect.
1625        let idx_res = energies
1626            .iter()
1627            .position(|&e| (e - 6.674).abs() < 0.05)
1628            .unwrap();
1629        assert!(
1630            dxs_dt[0][idx_res].abs() > 0.0,
1631            "dσ/dT should be non-zero near resonance, got {}",
1632            dxs_dt[0][idx_res]
1633        );
1634
1635        // Cross-check: compute a manual FD at a larger step and verify
1636        // the analytical derivative is consistent (within ~5% relative error).
1637        let big_dt = 1.0;
1638        let xs_up = broadened_cross_sections(
1639            &energies,
1640            std::slice::from_ref(&data),
1641            temperature + big_dt,
1642            None,
1643            None,
1644        )
1645        .unwrap();
1646        let xs_down =
1647            broadened_cross_sections(&energies, &[data], temperature - big_dt, None, None).unwrap();
1648
1649        let manual_deriv: Vec<f64> = xs_up[0]
1650            .iter()
1651            .zip(xs_down[0].iter())
1652            .map(|(&u, &d)| (u - d) / (2.0 * big_dt))
1653            .collect();
1654
1655        // Compare near the resonance where the derivative is large.
1656        let deriv_analytical = dxs_dt[0][idx_res];
1657        let deriv_coarse = manual_deriv[idx_res];
1658        let rel_err = (deriv_analytical - deriv_coarse).abs()
1659            / deriv_analytical.abs().max(deriv_coarse.abs()).max(1e-30);
1660        assert!(
1661            rel_err < 0.05,
1662            "Analytical vs FD derivatives disagree: analytical={}, coarse={}, rel_err={}",
1663            deriv_analytical,
1664            deriv_coarse,
1665            rel_err,
1666        );
1667    }
1668
1669    #[test]
1670    fn test_broadened_xs_analytical_derivative_low_temperature() {
1671        // Regression test: derivative must have correct sign at low temperature.
1672        let data = u238_single_resonance();
1673        let energies: Vec<f64> = (0..51).map(|i| 5.0 + (i as f64) * 0.1).collect();
1674
1675        let base_xs =
1676            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
1677
1678        // T = 0.05 K
1679        let (xs_low, dxs_low) = broadened_cross_sections_with_analytical_derivative_from_base(
1680            &energies,
1681            &base_xs,
1682            std::slice::from_ref(&data),
1683            0.05,
1684            None,
1685        )
1686        .unwrap();
1687        assert!(!xs_low.is_empty());
1688        // Derivative should be finite and mostly positive (Doppler broadening
1689        // increases with temperature for narrow resonances).
1690        for deriv_vec in &dxs_low {
1691            for &d in deriv_vec {
1692                assert!(d.is_finite(), "derivative must be finite at T=0.05 K");
1693            }
1694        }
1695
1696        // T = 0.0 K: edge case
1697        let (xs_zero, dxs_zero) = broadened_cross_sections_with_analytical_derivative_from_base(
1698            &energies,
1699            &base_xs,
1700            std::slice::from_ref(&data),
1701            0.0,
1702            None,
1703        )
1704        .unwrap();
1705        assert!(!xs_zero.is_empty());
1706        for deriv_vec in &dxs_zero {
1707            for &d in deriv_vec {
1708                assert!(d.is_finite(), "derivative must be finite at T=0.0 K");
1709            }
1710        }
1711    }
1712
1713    // --- SampleParams validation tests ---
1714
1715    #[test]
1716    fn test_sample_params_valid() {
1717        let sample = SampleParams::new(300.0, vec![]).unwrap();
1718        assert!((sample.temperature_k() - 300.0).abs() < 1e-15);
1719        assert!(sample.isotopes().is_empty());
1720    }
1721
1722    #[test]
1723    fn test_sample_params_zero_temperature() {
1724        let sample = SampleParams::new(0.0, vec![]).unwrap();
1725        assert!((sample.temperature_k()).abs() < 1e-15);
1726    }
1727
1728    #[test]
1729    fn test_sample_params_rejects_negative_temperature() {
1730        let err = SampleParams::new(-1.0, vec![]).unwrap_err();
1731        assert_eq!(err, SampleParamsError::NegativeTemperature(-1.0));
1732    }
1733
1734    #[test]
1735    fn test_sample_params_rejects_nan_temperature() {
1736        let err = SampleParams::new(f64::NAN, vec![]).unwrap_err();
1737        assert!(matches!(err, SampleParamsError::NonFiniteTemperature(_)));
1738    }
1739
1740    #[test]
1741    fn test_sample_params_rejects_infinite_temperature() {
1742        let err = SampleParams::new(f64::INFINITY, vec![]).unwrap_err();
1743        assert!(matches!(err, SampleParamsError::NonFiniteTemperature(_)));
1744    }
1745
1746    #[test]
1747    fn test_sample_params_rejects_neg_infinite_temperature() {
1748        let err = SampleParams::new(f64::NEG_INFINITY, vec![]).unwrap_err();
1749        assert!(matches!(err, SampleParamsError::NonFiniteTemperature(_)));
1750    }
1751
1752    // --- Base XS caching tests ---
1753
1754    #[test]
1755    fn test_forward_model_from_base_xs_matches_forward_model() {
1756        let data = u238_single_resonance();
1757        let thickness = 0.0005;
1758        let temperature = 300.0;
1759        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
1760
1761        // Reference: full forward model
1762        let sample = SampleParams::new(temperature, vec![(data.clone(), thickness)]).unwrap();
1763        let t_ref = forward_model(&energies, &sample, None).unwrap();
1764
1765        // Cached path: unbroadened XS → forward_model_from_base_xs
1766        let base_xs =
1767            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
1768        let t_cached = forward_model_from_base_xs(
1769            &energies,
1770            &base_xs,
1771            std::slice::from_ref(&data),
1772            &[thickness],
1773            temperature,
1774            None,
1775        )
1776        .unwrap();
1777
1778        for (i, (&r, &c)) in t_ref.iter().zip(t_cached.iter()).enumerate() {
1779            assert!(
1780                (r - c).abs() < 1e-12,
1781                "Mismatch at E[{}]={}: ref={}, cached={}",
1782                i,
1783                energies[i],
1784                r,
1785                c
1786            );
1787        }
1788    }
1789
1790    #[test]
1791    fn test_broadened_from_base_matches_broadened() {
1792        let data = u238_single_resonance();
1793        let temperature = 300.0;
1794        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
1795
1796        let xs_ref = broadened_cross_sections(
1797            &energies,
1798            std::slice::from_ref(&data),
1799            temperature,
1800            None,
1801            None,
1802        )
1803        .unwrap();
1804        let base_xs =
1805            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
1806        let xs_cached = broadened_cross_sections_from_base(
1807            &energies,
1808            &base_xs,
1809            std::slice::from_ref(&data),
1810            temperature,
1811            None,
1812        )
1813        .unwrap();
1814
1815        assert_eq!(xs_ref.len(), xs_cached.len());
1816        for (r, c) in xs_ref[0].iter().zip(xs_cached[0].iter()) {
1817            assert!(
1818                (r - c).abs() < 1e-12,
1819                "broadened_from_base mismatch: ref={}, cached={}",
1820                r,
1821                c
1822            );
1823        }
1824    }
1825
1826    #[test]
1827    fn test_analytical_derivative_from_base_shape_and_finiteness() {
1828        // Verify the analytical derivative from base returns correct shapes
1829        // and finite values.
1830        let data = u238_single_resonance();
1831        let temperature = 300.0;
1832        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
1833
1834        let base_xs =
1835            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
1836        let (xs, dxs_dt) = broadened_cross_sections_with_analytical_derivative_from_base(
1837            &energies,
1838            &base_xs,
1839            std::slice::from_ref(&data),
1840            temperature,
1841            None,
1842        )
1843        .unwrap();
1844
1845        assert_eq!(xs.len(), 1);
1846        assert_eq!(dxs_dt.len(), 1);
1847        assert_eq!(xs[0].len(), energies.len());
1848        assert_eq!(dxs_dt[0].len(), energies.len());
1849
1850        // All values should be finite.
1851        for &v in &xs[0] {
1852            assert!(v.is_finite(), "XS must be finite, got {v}");
1853        }
1854        for &v in &dxs_dt[0] {
1855            assert!(v.is_finite(), "dXS/dT must be finite, got {v}");
1856        }
1857
1858        // XS from base should match XS from full broadening.
1859        let xs_full = broadened_cross_sections(
1860            &energies,
1861            std::slice::from_ref(&data),
1862            temperature,
1863            None,
1864            None,
1865        )
1866        .unwrap();
1867        for (r, c) in xs_full[0].iter().zip(xs[0].iter()) {
1868            assert!(
1869                (r - c).abs() < 1e-12,
1870                "XS mismatch: full={}, from_base={}",
1871                r,
1872                c
1873            );
1874        }
1875    }
1876
1877    /// Regression test for issue #442: resolution broadening must be applied
1878    /// to the total transmission T(E) AFTER Beer-Lambert, not to σ(E) before.
1879    ///
1880    /// This test constructs the expected result from first principles:
1881    ///
1882    ///   1. Doppler-broaden σ
1883    ///   2. Beer-Lambert: T = exp(−n·σ_D)
1884    ///   3. Resolution-broaden T
1885    ///
1886    /// and asserts that `forward_model()` matches.
1887    #[test]
1888    fn test_forward_model_resolution_after_beer_lambert() {
1889        let data = u238_single_resonance();
1890        let thickness = 0.0005; // atoms/barn
1891        let temperature = 300.0;
1892
1893        // Energy grid around the 6.674 eV resonance.
1894        let energies: Vec<f64> = (0..401).map(|i| 4.0 + (i as f64) * 0.015).collect();
1895
1896        let inst = InstrumentParams {
1897            resolution: resolution::ResolutionFunction::Gaussian(
1898                resolution::ResolutionParams::new(25.0, 0.5, 0.005, 0.0).unwrap(),
1899            ),
1900        };
1901
1902        // --- Build expected from first principles ---
1903
1904        // Step 1: Doppler-broadened σ on the data grid.
1905        let unbroadened: Vec<f64> = energies
1906            .iter()
1907            .map(|&e| reich_moore::cross_sections_at_energy(&data, e).total)
1908            .collect();
1909        let doppler_params = doppler::DopplerParams::new(temperature, data.awr).unwrap();
1910        let sigma_d = doppler::doppler_broaden(&energies, &unbroadened, &doppler_params).unwrap();
1911
1912        // Step 2: Beer-Lambert on total transmission.
1913        let transmission: Vec<f64> = sigma_d.iter().map(|&s| (-thickness * s).exp()).collect();
1914
1915        // Step 3: Resolution-broaden the transmission.
1916        let t_expected =
1917            resolution::apply_resolution(&energies, &transmission, &inst.resolution).unwrap();
1918
1919        // --- Wrong ordering for comparison: Resolution(σ) then Beer-Lambert ---
1920        let sigma_broadened =
1921            resolution::apply_resolution(&energies, &sigma_d, &inst.resolution).unwrap();
1922        let t_wrong: Vec<f64> = sigma_broadened
1923            .iter()
1924            .map(|&s| (-thickness * s).exp())
1925            .collect();
1926
1927        // --- forward_model() output ---
1928        let sample = SampleParams::new(temperature, vec![(data, thickness)]).unwrap();
1929        let t_forward = forward_model(&energies, &sample, Some(&inst)).unwrap();
1930
1931        // forward_model should match the correct ordering (resolution after Beer-Lambert).
1932        // The extended grid in forward_model adds boundary points, so the match
1933        // is approximate — but should be very close on the interior grid.
1934        let interior = 20..energies.len() - 20; // skip boundary region
1935        let mut max_err_correct = 0.0f64;
1936        let mut max_err_wrong = 0.0f64;
1937        for i in interior.clone() {
1938            let err_correct = (t_forward[i] - t_expected[i]).abs();
1939            let err_wrong = (t_forward[i] - t_wrong[i]).abs();
1940            max_err_correct = max_err_correct.max(err_correct);
1941            max_err_wrong = max_err_wrong.max(err_wrong);
1942        }
1943
1944        // The key discriminant: forward_model must be much closer to the
1945        // correct ordering (resolution after Beer-Lambert) than to the wrong
1946        // ordering (resolution before Beer-Lambert).
1947        //
1948        // Small absolute differences (~1%) between forward_model and the
1949        // data-grid reference are expected because forward_model uses an
1950        // extended grid for boundary handling.
1951        assert!(
1952            max_err_correct < max_err_wrong,
1953            "forward_model is closer to the WRONG ordering than the correct one. \
1954             Error vs correct = {max_err_correct}, error vs wrong = {max_err_wrong}"
1955        );
1956
1957        // The error against the correct ordering should be at least 5× smaller
1958        // than the error against the wrong ordering.
1959        assert!(
1960            max_err_correct < max_err_wrong * 0.5,
1961            "forward_model should be clearly closer to the correct ordering. \
1962             Error vs correct = {max_err_correct}, error vs wrong = {max_err_wrong}, \
1963             ratio = {:.2}",
1964            max_err_correct / max_err_wrong
1965        );
1966
1967        // Verify the two orderings actually differ — if they don't, the test
1968        // is not exercising the bug.
1969        let ordering_diff: f64 = interior
1970            .map(|i| (t_expected[i] - t_wrong[i]).abs())
1971            .fold(0.0f64, f64::max);
1972        assert!(
1973            ordering_diff > 1e-4,
1974            "The two orderings should differ measurably at the resonance dip, \
1975             but max diff = {ordering_diff}. Test parameters may be too weak."
1976        );
1977    }
1978
1979    /// Issue #442 containment: `broadened_cross_sections()` must return
1980    /// Doppler-only σ even when `instrument` is `Some`.  Resolution
1981    /// broadening must NOT be applied inside this function.
1982    #[test]
1983    fn test_broadened_xs_is_doppler_only_with_instrument() {
1984        let data = u238_single_resonance();
1985        let temperature = 300.0;
1986        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
1987
1988        let inst = InstrumentParams {
1989            resolution: resolution::ResolutionFunction::Gaussian(
1990                resolution::ResolutionParams::new(25.0, 0.5, 0.005, 0.0).unwrap(),
1991            ),
1992        };
1993
1994        // With instrument (used for aux grid, but NOT for resolution broadening).
1995        let xs_with_inst = broadened_cross_sections(
1996            &energies,
1997            std::slice::from_ref(&data),
1998            temperature,
1999            Some(&inst),
2000            None,
2001        )
2002        .unwrap();
2003
2004        // Without instrument (pure Doppler on data grid).
2005        let xs_no_inst = broadened_cross_sections(
2006            &energies,
2007            std::slice::from_ref(&data),
2008            temperature,
2009            None,
2010            None,
2011        )
2012        .unwrap();
2013
2014        // Both should return Doppler-only σ.  The with-instrument path uses
2015        // the extended grid for Doppler which may produce slightly different
2016        // values, but they must be close (no resolution smoothing).
2017        assert_eq!(xs_with_inst.len(), 1);
2018        assert_eq!(xs_no_inst.len(), 1);
2019        assert_eq!(xs_with_inst[0].len(), energies.len());
2020
2021        // Compute what resolution-broadened σ would look like.
2022        let sigma_resolved =
2023            resolution::apply_resolution(&energies, &xs_no_inst[0], &inst.resolution).unwrap();
2024
2025        // The with-instrument result must NOT match the resolution-broadened version.
2026        // Near the resonance dip, resolution broadening smooths the peak — the
2027        // Doppler-only result should have a deeper dip than the resolved one.
2028        let idx_dip = energies
2029            .iter()
2030            .position(|&e| (e - 6.674).abs() < 0.05)
2031            .unwrap();
2032        let diff_doppler = (xs_with_inst[0][idx_dip] - xs_no_inst[0][idx_dip]).abs();
2033        let diff_resolved = (xs_with_inst[0][idx_dip] - sigma_resolved[idx_dip]).abs();
2034
2035        // The Doppler-only values from both paths should be closer to each
2036        // other than to the resolution-broadened value.
2037        assert!(
2038            diff_doppler < diff_resolved,
2039            "broadened_cross_sections with instrument should return Doppler-only σ, \
2040             not resolution-broadened σ.  \
2041             diff(with_inst, no_inst) = {diff_doppler}, \
2042             diff(with_inst, resolved) = {diff_resolved}"
2043        );
2044    }
2045
2046    /// Issue #442 containment: `broadened_cross_sections_from_base()` must
2047    /// return Doppler-only σ even when `instrument` is `Some`.
2048    #[test]
2049    fn test_broadened_from_base_is_doppler_only_with_instrument() {
2050        let data = u238_single_resonance();
2051        let temperature = 300.0;
2052        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
2053
2054        let inst = InstrumentParams {
2055            resolution: resolution::ResolutionFunction::Gaussian(
2056                resolution::ResolutionParams::new(25.0, 0.5, 0.005, 0.0).unwrap(),
2057            ),
2058        };
2059
2060        let base_xs =
2061            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
2062
2063        // With instrument.
2064        let xs_with_inst = broadened_cross_sections_from_base(
2065            &energies,
2066            &base_xs,
2067            std::slice::from_ref(&data),
2068            temperature,
2069            Some(&inst),
2070        )
2071        .unwrap();
2072
2073        // Without instrument.
2074        let xs_no_inst = broadened_cross_sections_from_base(
2075            &energies,
2076            &base_xs,
2077            std::slice::from_ref(&data),
2078            temperature,
2079            None,
2080        )
2081        .unwrap();
2082
2083        // Both should produce Doppler-only σ.  With-instrument may use
2084        // extended grid, but the extracted data-grid values should be close.
2085        let sigma_resolved =
2086            resolution::apply_resolution(&energies, &xs_no_inst[0], &inst.resolution).unwrap();
2087
2088        let idx_dip = energies
2089            .iter()
2090            .position(|&e| (e - 6.674).abs() < 0.05)
2091            .unwrap();
2092        let diff_doppler = (xs_with_inst[0][idx_dip] - xs_no_inst[0][idx_dip]).abs();
2093        let diff_resolved = (xs_with_inst[0][idx_dip] - sigma_resolved[idx_dip]).abs();
2094
2095        assert!(
2096            diff_doppler < diff_resolved,
2097            "broadened_cross_sections_from_base with instrument should return Doppler-only σ. \
2098             diff(with_inst, no_inst) = {diff_doppler}, \
2099             diff(with_inst, resolved) = {diff_resolved}"
2100        );
2101    }
2102
2103    // ── Issue #442 Step 5: forward_model_from_base_xs resolution ordering ──
2104
2105    /// Issue #442 Step 5: `forward_model_from_base_xs()` with resolution must
2106    /// match `forward_model()` with resolution for the same sample.
2107    #[test]
2108    fn test_forward_model_from_base_xs_matches_forward_model_with_resolution() {
2109        let data = u238_single_resonance();
2110        let thickness = 0.0005;
2111        let temperature = 300.0;
2112        let energies: Vec<f64> = (0..401).map(|i| 4.0 + (i as f64) * 0.015).collect();
2113
2114        let inst = InstrumentParams {
2115            resolution: resolution::ResolutionFunction::Gaussian(
2116                resolution::ResolutionParams::new(25.0, 0.5, 0.005, 0.0).unwrap(),
2117            ),
2118        };
2119
2120        // Reference: forward_model() (fixed in Step 1).
2121        let sample = SampleParams::new(temperature, vec![(data.clone(), thickness)]).unwrap();
2122        let t_ref = forward_model(&energies, &sample, Some(&inst)).unwrap();
2123
2124        // Base-XS path: unbroadened → forward_model_from_base_xs with resolution.
2125        let base_xs =
2126            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
2127        let t_base = forward_model_from_base_xs(
2128            &energies,
2129            &base_xs,
2130            std::slice::from_ref(&data),
2131            &[thickness],
2132            temperature,
2133            Some(&inst),
2134        )
2135        .unwrap();
2136
2137        // Both now run the IDENTICAL pipeline on the auxiliary extended grid
2138        // (Doppler → Beer-Lambert → resolution, data points extracted last), so
2139        // they must agree to floating-point round-off.  The cached path reuses
2140        // the data-grid base XS, which is bit-identical to forward_model's fresh
2141        // evaluation at those points; the auxiliary-only points are evaluated the
2142        // same way in both.  (Before Fix #7 the cached path applied resolution on
2143        // the coarse data grid and the two differed by ~1%.)
2144        let interior = 20..energies.len() - 20;
2145        let mut max_err = 0.0f64;
2146        for i in interior.clone() {
2147            max_err = max_err.max((t_ref[i] - t_base[i]).abs());
2148        }
2149        assert!(
2150            max_err < 1e-12,
2151            "forward_model_from_base_xs with resolution must match forward_model \
2152             to round-off (identical aux-grid pipeline).  Max error = {max_err}"
2153        );
2154
2155        // Verify resolution actually made a difference (not a vacuous test).
2156        let t_no_res = forward_model_from_base_xs(
2157            &energies,
2158            &base_xs,
2159            std::slice::from_ref(&data),
2160            &[thickness],
2161            temperature,
2162            None,
2163        )
2164        .unwrap();
2165        let res_diff: f64 = interior
2166            .map(|i| (t_base[i] - t_no_res[i]).abs())
2167            .fold(0.0f64, f64::max);
2168        assert!(
2169            res_diff > 1e-4,
2170            "Resolution should make a measurable difference, but max diff = {res_diff}"
2171        );
2172    }
2173
2174    /// Issue #442 Step 5: `forward_model_from_base_xs()` without resolution
2175    /// must remain unchanged (matches existing no-resolution test).
2176    #[test]
2177    fn test_forward_model_from_base_xs_no_resolution_unchanged() {
2178        let data = u238_single_resonance();
2179        let thickness = 0.0005;
2180        let temperature = 300.0;
2181        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
2182
2183        let sample = SampleParams::new(temperature, vec![(data.clone(), thickness)]).unwrap();
2184        let t_ref = forward_model(&energies, &sample, None).unwrap();
2185
2186        let base_xs =
2187            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
2188        let t_base = forward_model_from_base_xs(
2189            &energies,
2190            &base_xs,
2191            std::slice::from_ref(&data),
2192            &[thickness],
2193            temperature,
2194            None,
2195        )
2196        .unwrap();
2197
2198        for (i, (&r, &b)) in t_ref.iter().zip(t_base.iter()).enumerate() {
2199            assert!(
2200                (r - b).abs() < 1e-12,
2201                "No-resolution mismatch at E[{i}]={}: ref={r}, base={b}",
2202                energies[i]
2203            );
2204        }
2205    }
2206
2207    // ── Issue #442 Step 7: derivative helper containment ───────────────────
2208
2209    /// Issue #442 Step 7: `broadened_cross_sections_with_analytical_derivative_from_base()`
2210    /// must return Doppler-only σ and Doppler-only ∂σ/∂T even when instrument
2211    /// is present.  Resolution broadening must NOT be applied inside.
2212    #[test]
2213    fn test_derivative_helper_is_doppler_only_with_instrument() {
2214        let data = u238_single_resonance();
2215        let temperature = 300.0;
2216        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
2217
2218        let inst = InstrumentParams {
2219            resolution: resolution::ResolutionFunction::Gaussian(
2220                resolution::ResolutionParams::new(25.0, 0.5, 0.005, 0.0).unwrap(),
2221            ),
2222        };
2223
2224        let base_xs =
2225            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
2226
2227        // With instrument.
2228        let (xs_inst, dxs_inst) = broadened_cross_sections_with_analytical_derivative_from_base(
2229            &energies,
2230            &base_xs,
2231            std::slice::from_ref(&data),
2232            temperature,
2233            Some(&inst),
2234        )
2235        .unwrap();
2236
2237        // Without instrument.
2238        let (xs_none, dxs_none) = broadened_cross_sections_with_analytical_derivative_from_base(
2239            &energies,
2240            &base_xs,
2241            std::slice::from_ref(&data),
2242            temperature,
2243            None,
2244        )
2245        .unwrap();
2246
2247        assert_eq!(xs_inst.len(), 1);
2248        assert_eq!(dxs_inst.len(), 1);
2249
2250        // Both should return Doppler-only.  Compute what resolution-broadened
2251        // σ would look like, and verify the with-instrument result is closer
2252        // to the no-instrument result than to the resolved version.
2253        let sigma_resolved =
2254            resolution::apply_resolution(&energies, &xs_none[0], &inst.resolution).unwrap();
2255
2256        let idx_dip = energies
2257            .iter()
2258            .position(|&e| (e - 6.674).abs() < 0.05)
2259            .unwrap();
2260
2261        // σ check: with-instrument should be close to no-instrument (Doppler-only),
2262        // not to the resolution-broadened version.
2263        let diff_doppler = (xs_inst[0][idx_dip] - xs_none[0][idx_dip]).abs();
2264        let diff_resolved = (xs_inst[0][idx_dip] - sigma_resolved[idx_dip]).abs();
2265        assert!(
2266            diff_doppler < diff_resolved,
2267            "derivative helper σ with instrument should be Doppler-only. \
2268             diff(inst, none) = {diff_doppler}, diff(inst, resolved) = {diff_resolved}"
2269        );
2270
2271        // ∂σ/∂T check: same pattern — should be Doppler-only derivative.
2272        let dxs_resolved =
2273            resolution::apply_resolution(&energies, &dxs_none[0], &inst.resolution).unwrap();
2274        let ddiff_doppler = (dxs_inst[0][idx_dip] - dxs_none[0][idx_dip]).abs();
2275        let ddiff_resolved = (dxs_inst[0][idx_dip] - dxs_resolved[idx_dip]).abs();
2276        assert!(
2277            ddiff_doppler < ddiff_resolved,
2278            "derivative helper ∂σ/∂T with instrument should be Doppler-only. \
2279             diff(inst, none) = {ddiff_doppler}, diff(inst, resolved) = {ddiff_resolved}"
2280        );
2281    }
2282
2283    /// Issue #442 Step 7: derivative helper without resolution must be
2284    /// unchanged — same σ and ∂σ/∂T as before.
2285    #[test]
2286    fn test_derivative_helper_no_resolution_unchanged() {
2287        let data = u238_single_resonance();
2288        let temperature = 300.0;
2289        let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
2290
2291        let base_xs =
2292            unbroadened_cross_sections(&energies, std::slice::from_ref(&data), None).unwrap();
2293
2294        let (xs, dxs) = broadened_cross_sections_with_analytical_derivative_from_base(
2295            &energies,
2296            &base_xs,
2297            std::slice::from_ref(&data),
2298            temperature,
2299            None,
2300        )
2301        .unwrap();
2302
2303        // σ should match broadened_cross_sections (no instrument).
2304        let xs_ref = broadened_cross_sections(
2305            &energies,
2306            std::slice::from_ref(&data),
2307            temperature,
2308            None,
2309            None,
2310        )
2311        .unwrap();
2312
2313        for (i, (&a, &b)) in xs[0].iter().zip(xs_ref[0].iter()).enumerate() {
2314            assert!(
2315                (a - b).abs() < 1e-12,
2316                "σ mismatch at E[{i}]: derivative_helper={a}, broadened={b}"
2317            );
2318        }
2319
2320        // ∂σ/∂T should be finite and non-trivial near resonance.
2321        assert_eq!(dxs.len(), 1);
2322        assert_eq!(dxs[0].len(), energies.len());
2323        let idx_res = energies
2324            .iter()
2325            .position(|&e| (e - 6.674).abs() < 0.05)
2326            .unwrap();
2327        assert!(
2328            dxs[0][idx_res].abs() > 0.0,
2329            "∂σ/∂T should be non-zero near resonance"
2330        );
2331        for &d in &dxs[0] {
2332            assert!(d.is_finite(), "∂σ/∂T must be finite");
2333        }
2334    }
2335
2336    // ── Fix #8: broadened_cross_sections_for_transmission thickness guard ──
2337
2338    /// `broadened_cross_sections_for_transmission` divides by the thickness in
2339    /// the σ_eff = −ln(T)/nd inversion, so a non-positive or non-finite
2340    /// thickness must be rejected up-front (the docstring requires nd > 0).
2341    #[test]
2342    fn test_broadened_for_transmission_rejects_bad_thickness() {
2343        let data = u238_single_resonance();
2344        let energies: Vec<f64> = (0..51).map(|i| 5.0 + (i as f64) * 0.1).collect();
2345        let inst = InstrumentParams {
2346            resolution: resolution::ResolutionFunction::Gaussian(
2347                resolution::ResolutionParams::new(25.0, 0.5, 0.005, 0.0).unwrap(),
2348            ),
2349        };
2350
2351        for bad in [0.0, -1.0, f64::NAN, f64::INFINITY, f64::NEG_INFINITY] {
2352            let err = broadened_cross_sections_for_transmission(
2353                &energies,
2354                std::slice::from_ref(&data),
2355                300.0,
2356                &inst,
2357                bad,
2358                None,
2359            )
2360            .unwrap_err();
2361            assert!(
2362                matches!(err, TransmissionError::InputMismatch(_)),
2363                "thickness = {bad} should be rejected with InputMismatch, got {err:?}"
2364            );
2365        }
2366
2367        // A valid positive thickness still succeeds.
2368        let ok = broadened_cross_sections_for_transmission(
2369            &energies,
2370            std::slice::from_ref(&data),
2371            300.0,
2372            &inst,
2373            0.001,
2374            None,
2375        );
2376        assert!(ok.is_ok(), "valid thickness should succeed: {ok:?}");
2377    }
2378
2379    /// Issue #608 (PR #609 coverage): the working-grid `*_from_base` path's
2380    /// input validation — `validate_base_xs` shape errors and the unsorted-grid
2381    /// guard.  These error branches were untested (the temperature fit always
2382    /// supplies well-formed, sorted inputs).
2383    #[test]
2384    fn from_base_working_grid_rejects_malformed_inputs() {
2385        let rd = vec![u238_single_resonance()];
2386        let energies: Vec<f64> = (0..21).map(|i| 1.0 + (i as f64) * 0.1).collect();
2387        let n_e = energies.len();
2388        let good_base = vec![vec![10.0f64; n_e]];
2389        // base_xs isotope-count mismatch (2 rows vs 1 isotope).
2390        let e1 = broadened_cross_sections_with_analytical_derivative_from_base(
2391            &energies,
2392            &[vec![10.0; n_e], vec![10.0; n_e]],
2393            &rd,
2394            300.0,
2395            None,
2396        )
2397        .unwrap_err();
2398        assert!(e1.to_string().contains("isotopes"), "got: {e1}");
2399        // base_xs row-length mismatch.
2400        let e2 = broadened_cross_sections_with_analytical_derivative_from_base(
2401            &energies,
2402            &[vec![10.0; n_e - 1]],
2403            &rd,
2404            300.0,
2405            None,
2406        )
2407        .unwrap_err();
2408        assert!(e2.to_string().contains("energies"), "got: {e2}");
2409        // Unsorted energies + an instrument ⇒ the sorted-grid guard fires.
2410        let inst = InstrumentParams {
2411            resolution: resolution::ResolutionFunction::Gaussian(
2412                resolution::ResolutionParams::new(25.0, 0.5, 0.005, 0.0).unwrap(),
2413            ),
2414        };
2415        let mut unsorted = energies.clone();
2416        unsorted.swap(0, 1);
2417        let e3 = broadened_cross_sections_with_analytical_derivative_from_base(
2418            &unsorted,
2419            &good_base,
2420            &rd,
2421            300.0,
2422            Some(&inst),
2423        )
2424        .unwrap_err();
2425        assert!(e3.to_string().contains("sorted"), "got: {e3}");
2426    }
2427}