1use 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
37fn 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 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 n_pw_linear >= 3
79 } else {
80 true
81 };
82
83 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
105fn 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 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 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; }
137 if gd > 0.0 {
138 pairs.push((res.energy, gd));
139 }
140 }
141 }
142 }
143 }
144 }
145 pairs
146}
147
148pub 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 e.dedup();
171 e
172}
173
174fn 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
203fn 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
231fn 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
244fn 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
268pub 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#[derive(Debug)]
298pub enum TransmissionError {
299 Resolution(ResolutionError),
301 Doppler(DopplerParamsError),
303 DopplerBroadening(crate::doppler::DopplerError),
305 Cancelled,
307 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
353pub type BroadenedXsWithDerivative = (Vec<Vec<f64>>, Vec<Vec<f64>>);
359
360#[derive(Debug, Clone)]
382pub struct WorkingGridLayout {
383 pub energies: Vec<f64>,
386 pub data_indices: Vec<usize>,
389}
390
391impl WorkingGridLayout {
392 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 pub fn extract(&self, working: &[f64]) -> Vec<f64> {
405 self.data_indices.iter().map(|&i| working[i]).collect()
406 }
407}
408
409pub struct WorkingGridXs {
416 pub sigma: Vec<Vec<f64>>,
418 pub layout: WorkingGridLayout,
420}
421
422pub struct WorkingGridXsWithDerivative {
430 pub sigma: Vec<Vec<f64>>,
432 pub dsigma_dt: Vec<Vec<f64>>,
434 pub layout: WorkingGridLayout,
436}
437
438pub 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
455pub 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#[derive(Debug, PartialEq)]
509pub enum SampleParamsError {
510 NonFiniteTemperature(f64),
512 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#[derive(Debug, Clone)]
533pub struct SampleParams {
534 temperature_k: f64,
536 isotopes: Vec<(ResonanceData, f64)>,
538}
539
540impl SampleParams {
541 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 #[must_use]
565 pub fn temperature_k(&self) -> f64 {
566 self.temperature_k
567 }
568
569 #[must_use]
571 pub fn isotopes(&self) -> &[(ResonanceData, f64)] {
572 &self.isotopes
573 }
574}
575
576#[derive(Debug, Clone)]
578pub struct InstrumentParams {
579 pub resolution: ResolutionFunction,
581}
582
583pub 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 if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
620 return Err(ResolutionError::UnsortedEnergies.into());
621 }
622
623 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 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, ¶ms)?
673 } else {
674 unbroadened
675 };
676 Ok((after_doppler, *thickness))
677 })
678 .collect();
679 let doppler_xs = doppler_xs?;
680
681 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 let transmission: Vec<f64> = total_attenuation.iter().map(|&att| (-att).exp()).collect();
691
692 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
706pub 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 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
761pub 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 if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
779 return Err(ResolutionError::UnsortedEnergies.into());
780 }
781
782 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 let result: Result<Vec<Vec<f64>>, TransmissionError> = resonance_data
797 .par_iter()
798 .map(|rd| {
799 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, ¶ms).map_err(Into::into)
811 } else {
812 Ok(unbroadened)
813 }
814 })
815 .collect();
816
817 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
829pub 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 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 let unbroadened: Vec<f64> = ext_energies
895 .iter()
896 .map(|&e| reich_moore::cross_sections_at_energy(rd, e).total)
897 .collect();
898
899 let after_doppler = if temperature_k > 0.0 {
901 let params = DopplerParams::new(temperature_k, rd.awr)?;
902 doppler::doppler_broaden(ext_energies, &unbroadened, ¶ms)?
903 } else {
904 unbroadened
905 };
906
907 let transmission: Vec<f64> = after_doppler
909 .iter()
910 .map(|&sigma| (-nd * sigma).exp())
911 .collect();
912
913 let t_broadened = resolution::apply_resolution_presorted(
915 ext_energies,
916 &transmission,
917 &instrument.resolution,
918 );
919
920 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 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, ¶ms)?
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
974pub 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
1008pub 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 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
1039pub 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 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 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, ¶ms).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
1096pub 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 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
1132pub 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 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 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, ¶ms)
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 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
1198pub 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 if instrument.is_some() && !energies.windows(2).all(|w| w[0] <= w[1]) {
1254 return Err(ResolutionError::UnsortedEnergies.into());
1255 }
1256
1257 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 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 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, ¶ms).map_err(Into::into)
1291 } else {
1292 Ok(xs_ext)
1293 }
1294 })
1295 .collect();
1296 let doppler_xs = doppler_xs?;
1297
1298 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 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 #[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 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 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 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 let data = u238_single_resonance();
1432 let thickness = 0.001; 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 let t_on_res = trans[2];
1444 let t_off_res = trans[0]; 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 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 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 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 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 let data = u238_single_resonance();
1497 let thickness = 0.0001; let energies: Vec<f64> = (0..401).map(|i| 5.0 + (i as f64) * 0.01).collect();
1500
1501 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 let sample_hot = SampleParams::new(300.0, vec![(data, thickness)]).unwrap();
1507 let t_hot = forward_model(&energies, &sample_hot, None).unwrap();
1508
1509 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 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 let u238 = u238_single_resonance();
1526
1527 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 let idx_u238 = energies
1569 .iter()
1570 .position(|&e| (e - 6.7).abs() < 0.05)
1571 .unwrap();
1572 let idx_other = energies
1574 .iter()
1575 .position(|&e| (e - 20.0).abs() < 0.05)
1576 .unwrap();
1577 let idx_off = energies
1579 .iter()
1580 .position(|&e| (e - 15.0).abs() < 0.05)
1581 .unwrap();
1582
1583 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 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 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 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 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 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 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 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 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 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 #[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 #[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 let sample = SampleParams::new(temperature, vec![(data.clone(), thickness)]).unwrap();
1763 let t_ref = forward_model(&energies, &sample, None).unwrap();
1764
1765 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 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 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 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 #[test]
1888 fn test_forward_model_resolution_after_beer_lambert() {
1889 let data = u238_single_resonance();
1890 let thickness = 0.0005; let temperature = 300.0;
1892
1893 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 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 let transmission: Vec<f64> = sigma_d.iter().map(|&s| (-thickness * s).exp()).collect();
1914
1915 let t_expected =
1917 resolution::apply_resolution(&energies, &transmission, &inst.resolution).unwrap();
1918
1919 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 let sample = SampleParams::new(temperature, vec![(data, thickness)]).unwrap();
1929 let t_forward = forward_model(&energies, &sample, Some(&inst)).unwrap();
1930
1931 let interior = 20..energies.len() - 20; 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 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 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 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 #[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 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 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 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 let sigma_resolved =
2023 resolution::apply_resolution(&energies, &xs_no_inst[0], &inst.resolution).unwrap();
2024
2025 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 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 #[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 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 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 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 #[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 let sample = SampleParams::new(temperature, vec![(data.clone(), thickness)]).unwrap();
2122 let t_ref = forward_model(&energies, &sample, Some(&inst)).unwrap();
2123
2124 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 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 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 #[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 #[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 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 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 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 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 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 #[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 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 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 #[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 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 #[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 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 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 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}