nereids_physics/urr.rs
1//! Unresolved Resonance Region (LRU=2) cross-section calculation.
2//!
3//! Computes average cross-sections via the Hauser-Feshbach formula using
4//! the average level-spacing and partial widths stored in `UrrData`.
5//!
6//! ## Hauser-Feshbach Formula
7//!
8//! For each (L, J) combination at energy E:
9//!
10//! ```text
11//! g_J = (2J+1) / ((2I+1) · (2s+1))
12//!
13//! LRF=1: Γ_n(E) = 2 · P_L(ρ(E)) · GNO [GNO = reduced neutron width]
14//! LRF=2: Γ_n(E) from tabulated energy grid using INT interpolation
15//! INT=1 histogram, INT=2 lin-lin, INT=3 log-x/lin-y,
16//! INT=4 lin-x/log-y, INT=5 log-log (all 5 ENDF codes supported)
17//!
18//! Γ_tot = Γ_n + GG + GF + GX
19//!
20//! σ_γ += (π/k²) · g_J · (2π · Γ_n · GG) / (D · Γ_tot)
21//! σ_f += (π/k²) · g_J · (2π · Γ_n · GF) / (D · Γ_tot)
22//! σ_cn += (π/k²) · g_J · (2π · Γ_n · Γ_n) / (D · Γ_tot) [neutron-out]
23//! ```
24//!
25//! `σ_cn` is the compound-elastic contribution; the competitive (GX) channel
26//! contributes to `total` but not `elastic`. Potential scattering
27//! `σ_pot = 4π·AP²/100` (AP in fm, result in barns) is included in the
28//! returned `total` and `elastic` so that the URR band produces a physically
29//! consistent cross-section without requiring special handling at the call site.
30//!
31//! ## Width-Fluctuation Correction — not yet implemented
32//!
33//! The kernel above is the energy-averaged Hauser-Feshbach expression in the
34//! W = 1 limit: no width-fluctuation (Porter-Thomas / Moldauer / Dresner)
35//! factor is applied. The χ² degrees of freedom `AMUN` / `AMUF` are parsed
36//! from ENDF File 2 and stored on `UrrJGroup` (within `UrrData`), but are
37//! not yet consumed here.
38//! SAMMY's URR treatment (FITACS, manual Section VIII) includes the
39//! fluctuation integrals, so NEREIDS URR cross-sections will differ from
40//! SAMMY's where fluctuation effects are significant.
41//!
42//! ## Units
43//! All energies in eV, all lengths (AP, channel radii) in fm (true physics
44//! femtometers, 10⁻¹⁵ m), cross-sections in barns.
45//!
46//! ENDF stores radii in 10⁻¹² cm (= 10 fm); the parser converts to fm at
47//! parse time by multiplying by 10 (see `ENDF_RADIUS_TO_FM` in `parser.rs`),
48//! matching SAMMY's `FillSammyRmatrixFromRMat.cpp` line 422.
49//!
50//! ## SAMMY Reference
51//! - SAMMY manual Section VIII.A (equations for the Unresolved Resonance Region)
52//! - SAMMY's URR analysis code is FITACS (`acs/` in the SAMMY source tree),
53//! which includes the width-fluctuation corrections this module omits
54//!
55//! ## References
56//! - ENDF-6 Formats Manual §2.2.2
57//! - W. Hauser, H. Feshbach, Phys. Rev. 87 (1952) 366
58
59use nereids_endf::resonance::UrrData;
60
61use crate::channel;
62use crate::penetrability;
63
64/// Compute energy-averaged Hauser-Feshbach cross-sections in the Unresolved
65/// Resonance Region (width-fluctuation correction not yet applied; see the
66/// module docs).
67///
68/// Returns `(total, elastic, capture, fission)` in barns.
69/// All four are zero when `e_ev` falls outside the URR energy band
70/// `[urr.e_low, urr.e_high]`.
71///
72/// `elastic` = compound-nuclear elastic + smooth potential scattering.
73/// `total` = elastic + capture + fission + competitive (GX).
74/// Potential scattering `σ_pot = 4π·AP²/100` (AP in fm, result in barns)
75/// is folded in here rather than at the call site, so the URR band yields a
76/// physically complete elastic cross-section.
77///
78/// ## SAMMY Reference
79/// SAMMY manual Section VIII.A — energy-averaged cross-section equations.
80///
81/// # Arguments
82/// * `urr` — Parsed URR parameters (LRF=1 or LRF=2).
83/// * `e_ev` — Neutron lab-frame energy in eV.
84/// * `ap_fm` — Scattering radius in fm at this energy. The caller is
85/// responsible for evaluating any AP(E) table (NRO≠0) or falling back
86/// to the constant `urr.ap`.
87///
88/// # Panics
89/// Panics if `e_ev` is not finite or is non-positive. This is a
90/// defensive guard at the public boundary; the Python wrapper and the
91/// SAMMY-style dispatcher already validate the energy grid via
92/// `validate_energy_grid`, so this assertion only fires for direct
93/// callers (other Rust crates, tests) that bypass the grid check.
94pub fn urr_cross_sections(urr: &UrrData, e_ev: f64, ap_fm: f64) -> (f64, f64, f64, f64) {
95 // Defensive input validation at the public boundary (issue #558).
96 // The dominant call path through `cross_sections_at_energy` is already
97 // gated by ENDF range bounds and the Python wrapper's
98 // `validate_energy_grid`, so this assertion should never fire in
99 // production. It exists to make this `pub fn` safe for direct callers
100 // (other Rust crates, tests, future bindings) — surfacing malformed
101 // energies as a loud panic rather than letting NaN/∞ silently
102 // contaminate downstream arithmetic. One branch at the public entry
103 // is well outside any inner loop.
104 assert!(
105 e_ev.is_finite() && e_ev > 0.0,
106 "expected positive finite e_ev, got {e_ev}"
107 );
108
109 if e_ev < urr.e_low || e_ev > urr.e_high {
110 return (0.0, 0.0, 0.0, 0.0);
111 }
112
113 let mut sig_cap = 0.0_f64;
114 let mut sig_fiss = 0.0_f64;
115 let mut sig_compound_n = 0.0_f64;
116 let mut sig_competitive = 0.0_f64;
117
118 for lg in &urr.l_groups {
119 let awri = lg.awri;
120 let l = lg.l;
121
122 // Guard: non-positive AWRI makes k² = 0 → π/k² = inf.
123 // SAMMY skips URR entirely so never hits this; we must guard explicitly.
124 if awri <= 0.0 {
125 continue;
126 }
127
128 // π/k² in barns and channel parameter ρ = k·AP.
129 // Uses the same channel formulas as the RRR calculator for consistency.
130 let pi_over_k2 = channel::pi_over_k_squared_barns(e_ev, awri);
131 let rho = channel::rho(e_ev, awri, ap_fm);
132
133 // Blatt-Weisskopf penetrability P_L(ρ) for LRF=1 Γ_n calculation.
134 let p_l = penetrability::penetrability(l, rho);
135
136 for jg in &lg.j_groups {
137 // Statistical weight g_J = (2J+1) / ((2I+1)·(2s+1)).
138 // Target spin I = urr.spi, neutron spin s = 1/2.
139 let g_j = channel::statistical_weight(jg.j, urr.spi);
140
141 // Obtain effective widths at energy e_ev (module formula block).
142 //
143 // LRF=1: Γ_n(E) = 2·P_L(ρ(E))·GNO (GNO = reduced neutron width, eV);
144 // D, GG, GF are energy-independent (single-element vecs).
145 // LRF=2: widths interpolated from the energy table using the
146 // INT code stored per J-group; all five ENDF INT codes
147 // (1-5) are supported (see the dispatch below).
148 let (gn_eff, d_eff, gx_eff, gg_eff, gf_eff) = if urr.lrf == 1 {
149 let gno = jg.gn[0]; // reduced neutron width (eV)
150 // LFW=1: fission widths are energy-dependent (tabulated).
151 // LFW=0: fission width is a single constant.
152 let gf = if jg.gf.len() > 1 {
153 // Interpolate fission width from shared energy grid.
154 lin_lin_interp(&jg.energies, &jg.gf, e_ev)
155 } else {
156 jg.gf[0]
157 };
158 (2.0 * p_l * gno, jg.d[0], 0.0_f64, jg.gg[0], gf)
159 } else {
160 // LRF=2: dispatch on the stored interpolation law.
161 // ENDF-6 §0.5: INT=1 histogram, 2 lin-lin, 3 log-x/lin-y,
162 // 4 lin-x/log-y, 5 log-log.
163 let interp: fn(&[f64], &[f64], f64) -> f64 = match jg.int_code {
164 1 => histogram_interp,
165 3 => log_lin_interp,
166 4 => lin_log_interp,
167 5 => log_log_interp,
168 _ => lin_lin_interp, // INT=2 and any unknown → lin-lin
169 };
170 let gn_i = interp(&jg.energies, &jg.gn, e_ev);
171 let d_i = interp(&jg.energies, &jg.d, e_ev);
172 let gx_i = interp(&jg.energies, &jg.gx, e_ev);
173 let gg_i = interp(&jg.energies, &jg.gg, e_ev);
174 let gf_i = interp(&jg.energies, &jg.gf, e_ev);
175 (gn_i, d_i, gx_i, gg_i, gf_i)
176 };
177
178 // Skip degenerate entries.
179 if d_eff <= 0.0 || gn_eff <= 0.0 {
180 continue;
181 }
182
183 let gamma_tot = gn_eff + gg_eff + gf_eff + gx_eff;
184 if gamma_tot <= 0.0 {
185 continue;
186 }
187
188 // Hauser-Feshbach kernel:
189 // (π/k²) · g_J · (2π · Γ_n / D) · Γ_x / Γ_tot
190 //
191 // The 2π/D factor is the average level density times 2π
192 // (module formula block; manual Sec. VIII.A).
193 let prefactor = pi_over_k2 * g_j * (2.0 * std::f64::consts::PI * gn_eff) / d_eff;
194
195 sig_cap += prefactor * gg_eff / gamma_tot;
196 sig_fiss += prefactor * gf_eff / gamma_tot;
197 // Compound-elastic contribution: neutron absorbed into the compound
198 // nucleus and re-emitted as a neutron. The probability of re-emission
199 // into the neutron channel is Γ_n / Γ_tot, giving:
200 // σ_cn = (π/k²) · g_J · (2π·Γ_n/D) · Γ_n / Γ_tot
201 // Compound-elastic: Γ_n appears as both the entrance and the exit
202 // width (module formula block).
203 sig_compound_n += prefactor * gn_eff / gamma_tot;
204 // Competitive (inelastic) channel: GX is included in Γ_tot but goes
205 // into neither capture nor fission. It contributes to the total
206 // cross section without appearing in elastic, capture, or fission.
207 // GX enters Γ_tot only (module formula block).
208 sig_competitive += prefactor * gx_eff / gamma_tot;
209 }
210 }
211
212 // Smooth potential scattering: σ_pot = 4π·AP²
213 // AP is in fm; 1 barn = 100 fm².
214 // Folded into `elastic` per the function contract (see fn doc).
215 let sig_pot = 4.0 * std::f64::consts::PI * ap_fm * ap_fm / 100.0;
216 let elastic = sig_compound_n + sig_pot;
217 let total = elastic + sig_cap + sig_fiss + sig_competitive;
218 (total, elastic, sig_cap, sig_fiss)
219}
220
221/// Shared binary-search kernel for the two interpolation modes.
222///
223/// Handles empty-table, lower-clamp, and upper-clamp, then locates the
224/// bracketing interval and calls `interp_fn(x, x0, y0, x1, y1)`.
225/// `xs` must be strictly ascending and the same length as `ys`.
226#[inline]
227fn table_interp(
228 xs: &[f64],
229 ys: &[f64],
230 x: f64,
231 interp_fn: impl Fn(f64, f64, f64, f64, f64) -> f64,
232) -> f64 {
233 debug_assert_eq!(xs.len(), ys.len());
234 debug_assert!(
235 !xs.is_empty(),
236 "URR interpolation table must not be empty (parser bug?)"
237 );
238 if xs.is_empty() {
239 return 0.0;
240 }
241 if x <= xs[0] {
242 return ys[0];
243 }
244 if x >= xs[xs.len() - 1] {
245 return ys[xs.len() - 1];
246 }
247 let i = xs.partition_point(|&xi| xi <= x);
248 interp_fn(x, xs[i - 1], ys[i - 1], xs[i], ys[i])
249}
250
251/// Linear-linear interpolation (INT=2, clamped to table endpoints).
252fn lin_lin_interp(xs: &[f64], ys: &[f64], x: f64) -> f64 {
253 table_interp(xs, ys, x, |x, x0, y0, _x1, y1| {
254 let dx = _x1 - x0;
255 if dx.abs() < f64::EPSILON {
256 return y0;
257 }
258 y0 + (x - x0) / dx * (y1 - y0)
259 })
260}
261
262/// Histogram interpolation (INT=1, clamped to table endpoints).
263///
264/// Returns the y-value at the left endpoint of each interval (step function).
265/// ENDF-6 §0.5: "y = y_i for x_i ≤ x < x_{i+1}".
266fn histogram_interp(xs: &[f64], ys: &[f64], x: f64) -> f64 {
267 table_interp(xs, ys, x, |_x, _x0, y0, _x1, _y1| y0)
268}
269
270/// Log-x / Linear-y interpolation (INT=3, clamped to table endpoints).
271///
272/// y(x) = y₀ + [ln(x) - ln(x₀)] / [ln(x₁) - ln(x₀)] × (y₁ - y₀)
273/// Falls back to lin-lin when any x ≤ 0 (log undefined).
274fn log_lin_interp(xs: &[f64], ys: &[f64], x: f64) -> f64 {
275 table_interp(xs, ys, x, |x, x0, y0, x1, y1| {
276 if x <= 0.0 || x0 <= 0.0 || x1 <= 0.0 {
277 let dx = x1 - x0;
278 return if dx.abs() < f64::EPSILON {
279 y0
280 } else {
281 y0 + (x - x0) / dx * (y1 - y0)
282 };
283 }
284 let denom = x1.ln() - x0.ln();
285 if denom.abs() < f64::EPSILON {
286 return y0;
287 }
288 let t = (x.ln() - x0.ln()) / denom;
289 y0 + t * (y1 - y0)
290 })
291}
292
293/// Linear-x / Log-y interpolation (INT=4, clamped to table endpoints).
294///
295/// y(x) = y₀ × (y₁/y₀)^[(x - x₀) / (x₁ - x₀)]
296/// Falls back to lin-lin when any y ≤ 0 (log undefined).
297fn lin_log_interp(xs: &[f64], ys: &[f64], x: f64) -> f64 {
298 table_interp(xs, ys, x, |x, x0, y0, x1, y1| {
299 let dx = x1 - x0;
300 if dx.abs() < f64::EPSILON {
301 return y0;
302 }
303 if y0 <= 0.0 || y1 <= 0.0 {
304 return y0 + (x - x0) / dx * (y1 - y0);
305 }
306 let t = (x - x0) / dx;
307 (y0.ln() + t * (y1.ln() - y0.ln())).exp()
308 })
309}
310
311/// Log-log interpolation (INT=5, clamped to table endpoints).
312///
313/// Used for LRF=2 width tables in ENDF (common for heavy actinides such as
314/// U-238 where widths follow an approximate power-law energy dependence).
315/// Falls back to lin-lin when any bracket value is ≤ 0 (log undefined) or
316/// when the log-space denominator is too small for safe division.
317fn log_log_interp(xs: &[f64], ys: &[f64], x: f64) -> f64 {
318 table_interp(xs, ys, x, |x, x0, y0, x1, y1| {
319 let dx = x1 - x0;
320 // Guard against degenerate x-intervals (duplicate energy points).
321 if dx.abs() < f64::EPSILON {
322 return y0;
323 }
324 // Guard against non-positive values (log undefined).
325 if x <= 0.0 || x0 <= 0.0 || x1 <= 0.0 || y0 <= 0.0 || y1 <= 0.0 {
326 // Fall back to lin-lin for degenerate (non-positive) bracket entries.
327 return y0 + (x - x0) / dx * (y1 - y0);
328 }
329 let denom_ln = x1.ln() - x0.ln();
330 // Guard against nearly equal log values to avoid numerical blow-up.
331 if denom_ln.abs() < f64::EPSILON {
332 return y0;
333 }
334 let log_t = (x.ln() - x0.ln()) / denom_ln;
335 (y0.ln() + log_t * (y1.ln() - y0.ln())).exp()
336 })
337}
338
339#[cfg(test)]
340mod tests {
341 use super::*;
342 use nereids_endf::resonance::{UrrData, UrrJGroup, UrrLGroup};
343
344 /// σ = 0 when energy is outside the URR band.
345 #[test]
346 fn urr_outside_band_returns_zero() {
347 let urr = make_lrf1_urr(1000.0, 30_000.0);
348 let (t, e, c, f) = urr_cross_sections(&urr, 500.0, urr.ap); // below band
349 assert_eq!(t, 0.0);
350 assert_eq!(e, 0.0);
351 assert_eq!(c, 0.0);
352 assert_eq!(f, 0.0);
353
354 let (t2, _, _, _) = urr_cross_sections(&urr, 50_000.0, urr.ap); // above band
355 assert_eq!(t2, 0.0);
356 }
357
358 /// LRF=1 formula check: σ_γ = (π/k²) · g_J · (2π · Γ_n · GG) / (D · Γ_tot).
359 ///
360 /// Uses a single (L=0, J=2.0) group with known parameters; verifies
361 /// σ_γ matches the hand-computed value within numerical precision (~1e-10 relative).
362 #[test]
363 fn urr_lrf1_formula_check() {
364 // Hand-crafted parameters similar to U-233 URR (order of magnitude).
365 let e_test = 1_000.0_f64; // eV, mid-URR
366 let awri = 231.038_f64;
367 let ap = 9.6931_f64; // scattering radius in fm (ENDF 0.96931 × 10)
368 let spi = 2.5_f64; // U-233 target spin
369
370 let gno = 3.0e-4; // reduced neutron width (eV)
371 let gg = 3.5e-2; // gamma width (eV)
372 let d = 0.5; // level spacing (eV)
373 let j = 2.0_f64;
374
375 let urr = UrrData {
376 lrf: 1,
377 spi,
378 ap,
379 e_low: 100.0,
380 e_high: 30_000.0,
381 l_groups: vec![UrrLGroup {
382 l: 0,
383 awri,
384 j_groups: vec![UrrJGroup {
385 j,
386 amun: 1.0,
387 amuf: 0.0,
388 energies: vec![],
389 d: vec![d],
390 gx: vec![0.0],
391 gn: vec![gno],
392 gg: vec![gg],
393 gf: vec![0.0],
394 int_code: 2,
395 }],
396 }],
397 };
398
399 let (total, elastic, capture, fission) = urr_cross_sections(&urr, e_test, ap);
400
401 // Hand-compute expected values.
402 let pi_over_k2 = channel::pi_over_k_squared_barns(e_test, awri);
403 let rho = channel::rho(e_test, awri, ap);
404 let p_l = penetrability::penetrability(0, rho); // L=0
405 let gn_eff = 2.0 * p_l * gno;
406 let gamma_tot = gn_eff + gg;
407 let g_j = channel::statistical_weight(j, spi);
408 let prefactor = pi_over_k2 * g_j * (2.0 * std::f64::consts::PI * gn_eff) / d;
409 let expected_cap = prefactor * gg / gamma_tot;
410 let expected_cn = prefactor * gn_eff / gamma_tot; // neutron-out term
411 let expected_sig_pot = 4.0 * std::f64::consts::PI * ap * ap / 100.0;
412 let expected_elastic = expected_cn + expected_sig_pot;
413
414 assert!(
415 capture > 0.0,
416 "Capture cross-section must be positive, got {capture}"
417 );
418 assert_eq!(fission, 0.0, "No fission in this test case");
419 assert!(
420 (capture - expected_cap).abs() / expected_cap < 1e-10,
421 "Capture deviates from hand calculation: got {capture:.6e}, expected {expected_cap:.6e}"
422 );
423 assert!(
424 (elastic - expected_elastic).abs() / expected_elastic < 1e-10,
425 "Elastic deviates: got {elastic:.6e}, expected {expected_elastic:.6e} \
426 (compound_n={expected_cn:.6e} + sig_pot={expected_sig_pot:.6e})"
427 );
428 assert!(
429 (total - (capture + elastic)).abs() < 1e-14,
430 "Total ≠ capture + elastic: {total} ≠ {}",
431 capture + elastic
432 );
433 }
434
435 /// LRF=2 lin-lin interpolation smoke test.
436 ///
437 /// Two-point energy table; verifies σ_γ > 0 and changes monotonically
438 /// with energy (Γ_n increases with energy for L=0 neutrons).
439 #[test]
440 fn urr_lrf2_interpolation_smoke() {
441 let e1 = 1_000.0_f64;
442 let e2 = 10_000.0_f64;
443
444 let urr = UrrData {
445 lrf: 2,
446 spi: 2.5,
447 ap: 9.6931,
448 e_low: 600.0,
449 e_high: 30_000.0,
450 l_groups: vec![UrrLGroup {
451 l: 0,
452 awri: 231.038,
453 j_groups: vec![UrrJGroup {
454 j: 2.0,
455 amun: 1.0,
456 amuf: 4.0,
457 energies: vec![e1, e2],
458 d: vec![0.5, 0.4], // D decreases (more levels at higher E)
459 gx: vec![0.0, 0.0],
460 gn: vec![1.0e-3, 5.0e-3], // Γ_n increases with E
461 gg: vec![3.5e-2, 3.5e-2], // GG roughly constant
462 gf: vec![0.0, 0.0],
463 int_code: 2,
464 }],
465 }],
466 };
467
468 let (_, _, c1, _) = urr_cross_sections(&urr, e1, urr.ap);
469 let (_, _, c2, _) = urr_cross_sections(&urr, e2, urr.ap);
470 let (_, _, c_mid, _) = urr_cross_sections(&urr, 5_000.0, urr.ap);
471
472 assert!(c1 > 0.0, "σ_γ at {e1} eV must be positive, got {c1}");
473 assert!(c2 > 0.0, "σ_γ at {e2} eV must be positive, got {c2}");
474 // For a lin-lin table with two energy points, any interior energy
475 // should produce a σ_γ between the endpoint values.
476 let (c_lo, c_hi) = if c1 <= c2 { (c1, c2) } else { (c2, c1) };
477 assert!(
478 c_mid >= c_lo && c_mid <= c_hi,
479 "σ_γ at midpoint must lie between endpoint values; \
480 c1={c1:.3e} c_mid={c_mid:.3e} c2={c2:.3e}"
481 );
482
483 // Outside band returns zero.
484 let (t_out, _, _, _) = urr_cross_sections(&urr, 100.0, urr.ap);
485 assert_eq!(t_out, 0.0, "Outside band must be zero");
486 }
487
488 /// lin_lin_interp: boundary clamping and interior interpolation.
489 #[test]
490 fn lin_lin_interp_basic() {
491 let xs = vec![1.0, 3.0, 9.0];
492 let ys = vec![0.0, 2.0, 8.0];
493
494 // Below lower bound: clamp to y[0]
495 assert!((lin_lin_interp(&xs, &ys, 0.0) - 0.0).abs() < 1e-12);
496 // Above upper bound: clamp to y[last]
497 assert!((lin_lin_interp(&xs, &ys, 100.0) - 8.0).abs() < 1e-12);
498 // Interior: x=2.0, between [1,3], t=0.5 → y=1.0
499 assert!((lin_lin_interp(&xs, &ys, 2.0) - 1.0).abs() < 1e-12);
500 // Interior: x=6.0, between [3,9], t=0.5 → y=5.0
501 assert!((lin_lin_interp(&xs, &ys, 6.0) - 5.0).abs() < 1e-12);
502 }
503
504 /// log_log_interp: boundary clamping, power-law interior, and lin-lin fallback.
505 #[test]
506 fn log_log_interp_basic() {
507 // y = x^2: log-log should recover exact values.
508 let xs = vec![1.0, 10.0, 100.0];
509 let ys = vec![1.0, 100.0, 10_000.0]; // y = x^2
510
511 // Clamp below lower bound
512 assert!((log_log_interp(&xs, &ys, 0.5) - 1.0).abs() < 1e-12);
513 // Clamp above upper bound
514 assert!((log_log_interp(&xs, &ys, 200.0) - 10_000.0).abs() < 1e-12);
515 // Interior: x=sqrt(10)≈3.162, between [1,10]; exact power-law gives y=10
516 let x_mid = (1.0_f64 * 10.0_f64).sqrt(); // geometric mean → y = x_mid^2 = 10
517 assert!((log_log_interp(&xs, &ys, x_mid) - 10.0).abs() < 1e-10);
518 // At table points: exact values
519 assert!((log_log_interp(&xs, &ys, 10.0) - 100.0).abs() < 1e-10);
520 }
521
522 /// log_log_interp: falls back to lin-lin when a bracket value is ≤ 0.
523 #[test]
524 fn log_log_interp_nonpositive_fallback() {
525 // y[0] = 0 forces the lin-lin fallback in the first interval.
526 let xs = vec![1.0, 3.0, 9.0];
527 let ys = vec![0.0, 2.0, 8.0]; // y[0] = 0 → log undefined
528
529 // x=2 is in [1,3]; fallback: t=0.5, y = 0 + 0.5*(2-0) = 1.0
530 assert!((log_log_interp(&xs, &ys, 2.0) - 1.0).abs() < 1e-12);
531 // x=6 is in [3,9]; both y>0, so log-log applies (not a straight power law here,
532 // but result must be between 2 and 8)
533 let v = log_log_interp(&xs, &ys, 6.0);
534 assert!(
535 v > 2.0 && v < 8.0,
536 "log-log interior must be between endpoints: {v}"
537 );
538 }
539
540 // ─── Helpers ─────────────────────────────────────────────────────────────
541
542 fn make_lrf1_urr(e_low: f64, e_high: f64) -> UrrData {
543 UrrData {
544 lrf: 1,
545 spi: 2.5,
546 ap: 9.6931,
547 e_low,
548 e_high,
549 l_groups: vec![UrrLGroup {
550 l: 0,
551 awri: 231.038,
552 j_groups: vec![UrrJGroup {
553 j: 2.0,
554 amun: 1.0,
555 amuf: 0.0,
556 energies: vec![],
557 d: vec![0.5],
558 gx: vec![0.0],
559 gn: vec![1.0e-4],
560 gg: vec![3.5e-2],
561 gf: vec![0.0],
562 int_code: 2,
563 }],
564 }],
565 }
566 }
567
568 // ─── Defensive input validation at the public boundary ──────────────────
569 //
570 // `urr_cross_sections` is `pub`, so any Rust caller (other crates, tests,
571 // future integrations) can reach it directly without going through the
572 // Python wrapper's `validate_energy_grid`. The body's internal guards
573 // (`urr.e_low ≤ e_ev ≤ urr.e_high`) silently return zeros for negative
574 // energies because `e_low > 0` for any well-formed URR range — that
575 // would hide a caller bug. The entry assertion turns the footgun into
576 // a loud panic at the call site. See issue #558.
577
578 #[test]
579 #[should_panic(expected = "expected positive finite e_ev")]
580 fn urr_panics_on_zero_energy() {
581 let urr = make_lrf1_urr(1000.0, 30_000.0);
582 let _ = urr_cross_sections(&urr, 0.0, urr.ap);
583 }
584
585 #[test]
586 #[should_panic(expected = "expected positive finite e_ev")]
587 fn urr_panics_on_negative_energy() {
588 let urr = make_lrf1_urr(1000.0, 30_000.0);
589 let _ = urr_cross_sections(&urr, -1.0, urr.ap);
590 }
591
592 #[test]
593 #[should_panic(expected = "expected positive finite e_ev")]
594 fn urr_panics_on_nan_energy() {
595 let urr = make_lrf1_urr(1000.0, 30_000.0);
596 let _ = urr_cross_sections(&urr, f64::NAN, urr.ap);
597 }
598
599 #[test]
600 #[should_panic(expected = "expected positive finite e_ev")]
601 fn urr_panics_on_infinite_energy() {
602 let urr = make_lrf1_urr(1000.0, 30_000.0);
603 let _ = urr_cross_sections(&urr, f64::INFINITY, urr.ap);
604 }
605}