use std::sync::Arc;
use nereids_endf::resonance::ResonanceData;
use nereids_physics::transmission::{self, InstrumentParams};
use nereids_core::constants::{PIVOT_FLOOR, POISSON_EPSILON};
use crate::error::FittingError;
use crate::lm::FitModel;
use crate::parameters::ParameterSet;
#[derive(Debug, Clone)]
pub struct PoissonConfig {
pub max_iter: usize,
pub fd_step: f64,
pub step_size: f64,
pub tol_param: f64,
pub armijo_c: f64,
pub backtrack: f64,
pub lbfgsb_memory: usize,
}
impl Default for PoissonConfig {
fn default() -> Self {
Self {
max_iter: 200,
fd_step: 1e-7,
step_size: 1.0,
tol_param: 1e-8,
armijo_c: 1e-4,
backtrack: 0.5,
lbfgsb_memory: 7,
}
}
}
#[derive(Debug, Clone)]
pub struct PoissonResult {
pub nll: f64,
pub iterations: usize,
pub converged: bool,
pub params: Vec<f64>,
}
fn poisson_nll(y_obs: &[f64], y_model: &[f64]) -> f64 {
y_obs
.iter()
.zip(y_model.iter())
.map(|(&obs, &mdl)| poisson_nll_term(obs, mdl))
.sum()
}
#[inline]
fn poisson_nll_term(obs: f64, mdl: f64) -> f64 {
debug_assert!(
obs.is_finite() && obs >= 0.0,
"poisson_nll_term: obs must be finite and >= 0, got {obs}"
);
if mdl > POISSON_EPSILON {
mdl - obs * mdl.ln()
} else {
let eps = POISSON_EPSILON;
let nll_eps = eps - obs * eps.ln();
let grad_eps = 1.0 - obs / eps;
let hess_eps = if obs > 0.0 {
obs / (eps * eps)
} else {
1.0 / eps
};
let delta = eps - mdl;
nll_eps - grad_eps * delta + 0.5 * hess_eps * delta * delta
}
}
#[inline]
fn poisson_nll_grad_term(obs: f64, mdl: f64) -> f64 {
debug_assert!(
obs.is_finite() && obs >= 0.0,
"poisson_nll_grad_term: obs must be finite and >= 0, got {obs}"
);
if mdl > POISSON_EPSILON {
1.0 - obs / mdl
} else {
let eps = POISSON_EPSILON;
let grad_eps = 1.0 - obs / eps;
let hess_eps = if obs > 0.0 {
obs / (eps * eps)
} else {
1.0 / eps
};
let delta = eps - mdl;
grad_eps - hess_eps * delta
}
}
fn compute_gradient(
model: &dyn FitModel,
params: &mut ParameterSet,
y_obs: &[f64],
fd_step: f64,
all_vals_buf: &mut Vec<f64>,
free_idx_buf: &mut Vec<usize>,
) -> Result<Vec<f64>, FittingError> {
params.all_values_into(all_vals_buf);
let base_model = model.evaluate(all_vals_buf)?;
let base_nll = poisson_nll(y_obs, &base_model);
params.free_indices_into(free_idx_buf);
let mut grad = vec![0.0; free_idx_buf.len()];
for (j, &idx) in free_idx_buf.iter().enumerate() {
let original = params.params[idx].value;
let step = fd_step * (1.0 + original.abs());
params.params[idx].value = original + step;
params.params[idx].clamp();
let mut actual_step = params.params[idx].value - original;
if actual_step.abs() < PIVOT_FLOOR {
params.params[idx].value = original - step;
params.params[idx].clamp();
actual_step = params.params[idx].value - original;
if actual_step.abs() < PIVOT_FLOOR {
params.params[idx].value = original;
continue;
}
}
params.all_values_into(all_vals_buf);
let perturbed_model = match model.evaluate(all_vals_buf) {
Ok(v) => v,
Err(_) => {
params.params[idx].value = original;
continue;
}
};
let perturbed_nll = poisson_nll(y_obs, &perturbed_model);
params.params[idx].value = original;
grad[j] = (perturbed_nll - base_nll) / actual_step;
}
Ok(grad)
}
fn project(params: &mut ParameterSet) {
for p in &mut params.params {
p.clamp();
}
}
#[allow(clippy::too_many_arguments)]
fn backtracking_line_search(
model: &dyn FitModel,
params: &mut ParameterSet,
y_obs: &[f64],
old_free: &[f64],
search_dir: &[f64],
initial_alpha: f64,
config: &PoissonConfig,
grad: &[f64],
nll: f64,
all_vals_buf: &mut Vec<f64>,
free_vals_buf: &mut Vec<f64>,
trial_free_buf: &mut Vec<f64>,
) -> Option<f64> {
let mut alpha = initial_alpha;
for _ in 0..50 {
trial_free_buf.clear();
trial_free_buf.extend(
old_free
.iter()
.zip(search_dir.iter())
.map(|(&v, &d)| v - alpha * d),
);
params.set_free_values(trial_free_buf);
project(params);
params.all_values_into(all_vals_buf);
let trial_model = match model.evaluate(all_vals_buf) {
Ok(v) => v,
Err(_) => {
alpha *= config.backtrack;
continue;
}
};
if trial_model.iter().any(|v| !v.is_finite()) {
alpha *= config.backtrack;
continue;
}
let trial_nll = poisson_nll(y_obs, &trial_model);
params.free_values_into(free_vals_buf);
let descent = grad
.iter()
.zip(old_free.iter())
.zip(free_vals_buf.iter())
.map(|((&g, &old), &new)| g * (old - new))
.sum::<f64>();
if trial_nll.is_finite() && trial_nll <= nll - config.armijo_c * descent {
return Some(trial_nll);
}
alpha *= config.backtrack;
}
params.set_free_values(old_free);
None
}
fn try_early_return_fixed(
model: &dyn FitModel,
y_obs: &[f64],
params: &ParameterSet,
) -> Result<Option<PoissonResult>, FittingError> {
if params.n_free() != 0 {
return Ok(None);
}
let y_model = model.evaluate(¶ms.all_values())?;
let nll = poisson_nll(y_obs, &y_model);
if !nll.is_finite() {
return Ok(Some(PoissonResult {
nll,
iterations: 0,
converged: false,
params: params.all_values(),
}));
}
Ok(Some(PoissonResult {
nll,
iterations: 0,
converged: true,
params: params.all_values(),
}))
}
fn validate_analytic_inputs(
n_e: usize,
flux: &[f64],
cross_sections: &[Vec<f64>],
density_indices: &[usize],
params: &ParameterSet,
temp_ctx: Option<&TemperatureContext>,
caller: &str,
) -> Result<(), FittingError> {
if flux.len() != n_e {
return Err(FittingError::LengthMismatch {
expected: n_e,
actual: flux.len(),
field: "flux",
});
}
if density_indices.len() != cross_sections.len() {
return Err(FittingError::LengthMismatch {
expected: cross_sections.len(),
actual: density_indices.len(),
field: "density_indices",
});
}
for (k, sigma) in cross_sections.iter().enumerate() {
if sigma.len() != n_e {
return Err(FittingError::LengthMismatch {
expected: n_e,
actual: sigma.len(),
field: if k == 0 {
"cross_sections[0]"
} else {
"cross_sections[k]"
},
});
}
}
let free_set: std::collections::HashSet<usize> = params.free_indices().into_iter().collect();
let mut mapped_set: std::collections::HashSet<usize> =
density_indices.iter().copied().collect();
if let Some(ctx) = temp_ctx {
mapped_set.insert(ctx.temperature_index);
}
let unmapped: Vec<usize> = free_set.difference(&mapped_set).copied().collect();
if !unmapped.is_empty() {
return Err(FittingError::InvalidConfig(format!(
"{caller}: free parameters {unmapped:?} are not mapped by density_indices \
or temperature_index; analytical gradient is zero for these params. \
Use poisson_fit (FD) for mixed parameter sets."
)));
}
Ok(())
}
fn build_param_isotope_map(free_indices: &[usize], density_indices: &[usize]) -> Vec<Vec<usize>> {
free_indices
.iter()
.map(|&pi| {
density_indices
.iter()
.enumerate()
.filter(|&(_, &di)| di == pi)
.map(|(k, _)| k)
.collect()
})
.collect()
}
fn compute_transmission(
n_e: usize,
xs: &[Vec<f64>],
density_indices: &[usize],
all_vals: &[f64],
) -> Vec<f64> {
debug_assert!(
density_indices.iter().all(|&idx| idx < all_vals.len()),
"density_indices out of bounds for parameter vector"
);
(0..n_e)
.map(|e| {
let mut neg_opt = 0.0f64;
for (k, xs_k) in xs.iter().enumerate() {
let density = all_vals[density_indices[k]];
neg_opt -= density * xs_k[e];
}
neg_opt.exp()
})
.collect()
}
struct AnalyticGradientCtx<'a> {
y_obs: &'a [f64],
y_model: &'a [f64],
flux: &'a [f64],
t_now: &'a [f64],
xs: &'a [Vec<f64>],
param_isotopes: &'a [Vec<usize>],
temp_free_pos: Option<usize>,
density_indices: &'a [usize],
all_vals: &'a [f64],
dxs_dt: &'a [Vec<f64>],
}
fn compute_analytic_gradient(ctx: &AnalyticGradientCtx<'_>) -> Vec<f64> {
let mut grad: Vec<f64> = ctx
.param_isotopes
.iter()
.map(|iso_indices| {
ctx.y_obs
.iter()
.zip(ctx.y_model.iter())
.zip(ctx.flux.iter())
.zip(ctx.t_now.iter())
.enumerate()
.map(|(e, (((&obs, &ym), &phi), &t))| {
let residual_factor = poisson_nll_grad_term(obs, ym);
let sigma_sum: f64 = iso_indices.iter().map(|&k| ctx.xs[k][e]).sum();
residual_factor * phi * t * (-sigma_sum)
})
.sum::<f64>()
})
.collect();
if let Some(pos) = ctx.temp_free_pos {
debug_assert!(
!ctx.dxs_dt.is_empty(),
"dxs_dt must be non-empty when temperature fitting is enabled"
);
let temp_grad: f64 = ctx
.y_obs
.iter()
.zip(ctx.y_model.iter())
.zip(ctx.flux.iter())
.zip(ctx.t_now.iter())
.enumerate()
.map(|(e, (((&obs, &ym), &phi), &t))| {
let residual_factor = poisson_nll_grad_term(obs, ym);
let dsigma_sum: f64 = (0..ctx.density_indices.len())
.map(|k| {
let density = ctx.all_vals[ctx.density_indices[k]];
density * ctx.dxs_dt[k][e]
})
.sum();
residual_factor * phi * t * (-dsigma_sum)
})
.sum();
grad[pos] = temp_grad;
}
grad
}
pub fn poisson_fit(
model: &dyn FitModel,
y_obs: &[f64],
params: &mut ParameterSet,
config: &PoissonConfig,
) -> Result<PoissonResult, FittingError> {
if let Some(result) = try_early_return_fixed(model, y_obs, params)? {
return Ok(result);
}
let mut all_vals_buf = Vec::with_capacity(params.params.len());
let mut free_vals_buf = Vec::with_capacity(params.n_free());
let mut old_free_buf: Vec<f64> = Vec::with_capacity(params.n_free());
let mut trial_free_buf: Vec<f64> = Vec::with_capacity(params.n_free());
let mut free_idx_buf: Vec<usize> = Vec::with_capacity(params.n_free());
params.all_values_into(&mut all_vals_buf);
let y_model = model.evaluate(&all_vals_buf)?;
let mut nll = poisson_nll(y_obs, &y_model);
if !nll.is_finite() {
return Ok(PoissonResult {
nll,
iterations: 0,
converged: false,
params: params.all_values(),
});
}
let mut converged = false;
let mut iter = 0;
for _ in 0..config.max_iter {
iter += 1;
let grad = compute_gradient(
model,
params,
y_obs,
config.fd_step,
&mut all_vals_buf,
&mut free_idx_buf,
)?;
let grad_norm: f64 = grad.iter().map(|g| g * g).sum::<f64>().sqrt();
if grad_norm < config.tol_param {
converged = true;
break;
}
params.free_values_into(&mut free_vals_buf);
old_free_buf.clear();
old_free_buf.extend_from_slice(&free_vals_buf);
let initial_alpha = config.step_size / grad_norm.max(1.0);
match backtracking_line_search(
model,
params,
y_obs,
&old_free_buf,
&grad,
initial_alpha,
config,
&grad,
nll,
&mut all_vals_buf,
&mut free_vals_buf,
&mut trial_free_buf,
) {
Some(new_nll) => nll = new_nll,
None => {
break;
}
}
params.free_values_into(&mut free_vals_buf);
let step_norm: f64 = old_free_buf
.iter()
.zip(free_vals_buf.iter())
.map(|(o, n)| (o - n).powi(2))
.sum::<f64>()
.sqrt();
if step_norm < config.tol_param {
converged = true;
break;
}
}
Ok(PoissonResult {
nll,
iterations: iter,
converged,
params: params.all_values(),
})
}
pub struct TemperatureContext {
pub temperature_index: usize,
pub resonance_data: Vec<ResonanceData>,
pub energies: Vec<f64>,
pub instrument: Option<InstrumentParams>,
pub base_xs: Option<Arc<Vec<Vec<f64>>>>,
}
#[allow(clippy::too_many_arguments)]
pub fn poisson_fit_analytic(
model: &dyn FitModel,
y_obs: &[f64],
flux: &[f64],
cross_sections: &[Vec<f64>],
density_indices: &[usize],
params: &mut ParameterSet,
config: &PoissonConfig,
temp_ctx: Option<&TemperatureContext>,
) -> Result<PoissonResult, FittingError> {
let n_e = y_obs.len();
validate_analytic_inputs(
n_e,
flux,
cross_sections,
density_indices,
params,
temp_ctx,
"poisson_fit_analytic",
)?;
if let Some(result) = try_early_return_fixed(model, y_obs, params)? {
return Ok(result);
}
let free_indices = params.free_indices();
let param_isotopes = build_param_isotope_map(&free_indices, density_indices);
let temp_free_pos: Option<usize> = temp_ctx.as_ref().and_then(|ctx| {
free_indices
.iter()
.position(|&fi| fi == ctx.temperature_index)
});
let mut xs_owned: Vec<Vec<f64>> = if temp_ctx.is_some() {
cross_sections.to_vec()
} else {
Vec::new()
};
let mut dxs_dt: Vec<Vec<f64>> = if temp_ctx.is_some() {
vec![vec![0.0; n_e]; density_indices.len()]
} else {
Vec::new()
};
let mut all_vals_buf = Vec::with_capacity(params.params.len());
let mut free_vals_buf = Vec::with_capacity(params.n_free());
let mut old_free_buf: Vec<f64> = Vec::with_capacity(params.n_free());
let mut trial_free_buf: Vec<f64> = Vec::with_capacity(params.n_free());
params.all_values_into(&mut all_vals_buf);
let y_model = model.evaluate(&all_vals_buf)?;
let mut nll = poisson_nll(y_obs, &y_model);
if !nll.is_finite() {
return Ok(PoissonResult {
nll,
iterations: 0,
converged: false,
params: params.all_values(),
});
}
let mut converged = false;
let mut iter = 0;
let mut last_temperature: f64 = f64::NAN;
for _ in 0..config.max_iter {
iter += 1;
if let Some(ctx) = &temp_ctx {
params.all_values_into(&mut all_vals_buf);
let t_current = all_vals_buf[ctx.temperature_index];
if !last_temperature.is_finite() || (t_current - last_temperature).abs() > 1e-15 {
let (xs_new, dxs_new) = if let Some(ref base) = ctx.base_xs {
transmission::broadened_cross_sections_with_derivative_from_base(
&ctx.energies,
base,
&ctx.resonance_data,
t_current,
ctx.instrument.as_ref(),
)
} else {
transmission::broadened_cross_sections_with_derivative(
&ctx.energies,
&ctx.resonance_data,
t_current,
ctx.instrument.as_ref(),
)
}
.map_err(|e| {
FittingError::EvaluationFailed(format!(
"poisson_fit_analytic: broadening failed: {e}"
))
})?;
xs_owned = xs_new;
dxs_dt = dxs_new;
last_temperature = t_current;
}
}
let xs_ref: &[Vec<f64>] = if temp_ctx.is_some() {
&xs_owned
} else {
cross_sections
};
params.all_values_into(&mut all_vals_buf);
let y_model_now = model.evaluate(&all_vals_buf)?;
let t_now = compute_transmission(n_e, xs_ref, density_indices, &all_vals_buf);
let grad = compute_analytic_gradient(&AnalyticGradientCtx {
y_obs,
y_model: &y_model_now,
flux,
t_now: &t_now,
xs: xs_ref,
param_isotopes: ¶m_isotopes,
temp_free_pos,
density_indices,
all_vals: &all_vals_buf,
dxs_dt: &dxs_dt,
});
let grad_norm: f64 = grad.iter().map(|g| g * g).sum::<f64>().sqrt();
if grad_norm < config.tol_param {
converged = true;
break;
}
params.free_values_into(&mut free_vals_buf);
old_free_buf.clear();
old_free_buf.extend_from_slice(&free_vals_buf);
let hessian_diag: Vec<f64> = param_isotopes
.iter()
.enumerate()
.map(|(j, iso_indices)| {
let is_temp = temp_free_pos == Some(j);
y_model_now
.iter()
.zip(flux.iter())
.zip(t_now.iter())
.enumerate()
.map(|(e, ((&ym, &phi), &t))| {
let ym_safe = ym.max(POISSON_EPSILON);
let dy = if is_temp {
let dsigma_sum: f64 = (0..density_indices.len())
.map(|k| {
let density = all_vals_buf[density_indices[k]];
density * dxs_dt[k][e]
})
.sum();
phi * t * (-dsigma_sum)
} else {
let sigma_sum: f64 = iso_indices.iter().map(|&k| xs_ref[k][e]).sum();
phi * t * (-sigma_sum)
};
dy * dy / ym_safe
})
.sum::<f64>()
})
.collect();
let search_dir: Vec<f64> = grad
.iter()
.zip(hessian_diag.iter())
.map(|(&g, &h)| {
if h > PIVOT_FLOOR {
g / h
} else {
g
}
})
.collect();
let search_norm: f64 = search_dir.iter().map(|g| g * g).sum::<f64>().sqrt();
let initial_alpha = config.step_size / search_norm.max(1.0);
match backtracking_line_search(
model,
params,
y_obs,
&old_free_buf,
&search_dir,
initial_alpha,
config,
&grad,
nll,
&mut all_vals_buf,
&mut free_vals_buf,
&mut trial_free_buf,
) {
Some(new_nll) => nll = new_nll,
None => {
break;
}
}
params.free_values_into(&mut free_vals_buf);
let step_norm: f64 = old_free_buf
.iter()
.zip(free_vals_buf.iter())
.map(|(o, n)| (o - n).powi(2))
.sum::<f64>()
.sqrt();
if step_norm < config.tol_param {
converged = true;
break;
}
}
Ok(PoissonResult {
nll,
iterations: iter,
converged,
params: params.all_values(),
})
}
pub struct CountsModel<'a> {
pub transmission_model: &'a dyn FitModel,
pub flux: &'a [f64],
pub background: &'a [f64],
}
impl<'a> FitModel for CountsModel<'a> {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let transmission = self.transmission_model.evaluate(params)?;
debug_assert_eq!(
transmission.len(),
self.flux.len(),
"CountsModel: transmission length ({}) != flux length ({})",
transmission.len(),
self.flux.len(),
);
debug_assert_eq!(
self.flux.len(),
self.background.len(),
"CountsModel: flux length ({}) != background length ({})",
self.flux.len(),
self.background.len(),
);
Ok(transmission
.iter()
.zip(self.flux.iter())
.zip(self.background.iter())
.map(|((&t, &f), &b)| f * t + b)
.collect())
}
}
struct LbfgsbMemory {
m: usize,
s: Vec<Vec<f64>>,
y: Vec<Vec<f64>>,
rho: Vec<f64>,
len: usize,
cursor: usize,
}
impl LbfgsbMemory {
fn new(m: usize, n: usize) -> Self {
Self {
m,
s: vec![vec![0.0; n]; m],
y: vec![vec![0.0; n]; m],
rho: vec![0.0; m],
len: 0,
cursor: 0,
}
}
fn push(&mut self, s_vec: Vec<f64>, y_vec: Vec<f64>) {
let ys: f64 = y_vec
.iter()
.zip(s_vec.iter())
.map(|(&yi, &si)| yi * si)
.sum();
if ys <= PIVOT_FLOOR {
return;
}
let idx = self.cursor;
self.s[idx] = s_vec;
self.y[idx] = y_vec;
self.rho[idx] = 1.0 / ys;
self.cursor = (self.cursor + 1) % self.m;
if self.len < self.m {
self.len += 1;
}
}
fn apply_two_loop(&self, q: &[f64]) -> Vec<f64> {
let mut r: Vec<f64> = q.to_vec();
if self.len == 0 {
return r;
}
let mut alpha = vec![0.0; self.len];
for i in (0..self.len).rev() {
let idx = (self.cursor + self.m - self.len + i) % self.m;
let dot_sr: f64 = self.s[idx]
.iter()
.zip(r.iter())
.map(|(&si, &ri)| si * ri)
.sum();
alpha[i] = self.rho[idx] * dot_sr;
for (rj, &yj) in r.iter_mut().zip(self.y[idx].iter()) {
*rj -= alpha[i] * yj;
}
}
let newest = (self.cursor + self.m - 1) % self.m;
let sy: f64 = self.s[newest]
.iter()
.zip(self.y[newest].iter())
.map(|(&si, &yi)| si * yi)
.sum();
let yy: f64 = self.y[newest].iter().map(|&yi| yi * yi).sum();
if yy > PIVOT_FLOOR {
let gamma = sy / yy;
for ri in r.iter_mut() {
*ri *= gamma;
}
}
for (i, &alpha_i) in alpha.iter().enumerate().take(self.len) {
let idx = (self.cursor + self.m - self.len + i) % self.m;
let dot_yr: f64 = self.y[idx]
.iter()
.zip(r.iter())
.map(|(&yi, &ri)| yi * ri)
.sum();
let beta = self.rho[idx] * dot_yr;
for (rj, &sj) in r.iter_mut().zip(self.s[idx].iter()) {
*rj += (alpha_i - beta) * sj;
}
}
r
}
}
fn is_descent_direction(direction: &[f64], gradient: &[f64]) -> bool {
let dot: f64 = direction
.iter()
.zip(gradient.iter())
.map(|(&d, &g)| d * g)
.sum();
dot > 0.0
}
#[allow(clippy::too_many_arguments)]
pub fn poisson_fit_lbfgsb(
model: &dyn FitModel,
y_obs: &[f64],
flux: &[f64],
cross_sections: &[Vec<f64>],
density_indices: &[usize],
params: &mut ParameterSet,
config: &PoissonConfig,
temp_ctx: Option<&TemperatureContext>,
) -> Result<PoissonResult, FittingError> {
let n_e = y_obs.len();
if config.lbfgsb_memory < 1 {
return Err(FittingError::InvalidConfig(format!(
"lbfgsb_memory must be >= 1, got {}",
config.lbfgsb_memory,
)));
}
validate_analytic_inputs(
n_e,
flux,
cross_sections,
density_indices,
params,
temp_ctx,
"poisson_fit_lbfgsb",
)?;
if let Some(result) = try_early_return_fixed(model, y_obs, params)? {
return Ok(result);
}
let free_indices = params.free_indices();
let n_free = free_indices.len();
let param_isotopes = build_param_isotope_map(&free_indices, density_indices);
let temp_free_pos: Option<usize> = temp_ctx.as_ref().and_then(|ctx| {
free_indices
.iter()
.position(|&fi| fi == ctx.temperature_index)
});
let mut xs_owned: Vec<Vec<f64>> = if temp_ctx.is_some() {
cross_sections.to_vec()
} else {
Vec::new()
};
let mut dxs_dt: Vec<Vec<f64>> = if temp_ctx.is_some() {
vec![vec![0.0; n_e]; density_indices.len()]
} else {
Vec::new()
};
let mut all_vals_buf = Vec::with_capacity(params.params.len());
let mut free_vals_buf = Vec::with_capacity(params.n_free());
let mut trial_free_buf: Vec<f64> = Vec::with_capacity(params.n_free());
params.all_values_into(&mut all_vals_buf);
let y_model = model.evaluate(&all_vals_buf)?;
let mut nll = poisson_nll(y_obs, &y_model);
if !nll.is_finite() {
return Ok(PoissonResult {
nll,
iterations: 0,
converged: false,
params: params.all_values(),
});
}
let mut converged = false;
let mut iter = 0;
let mut last_temperature: f64 = f64::NAN;
let mut lbfgs = LbfgsbMemory::new(config.lbfgsb_memory, n_free);
params.free_values_into(&mut free_vals_buf);
let mut prev_free: Vec<f64> = free_vals_buf.clone();
let mut prev_grad: Option<Vec<f64>> = None;
for _ in 0..config.max_iter {
iter += 1;
params.all_values_into(&mut all_vals_buf);
if let Some(ctx) = &temp_ctx {
let t_current = all_vals_buf[ctx.temperature_index];
if !last_temperature.is_finite() || (t_current - last_temperature).abs() > 1e-15 {
let (xs_new, dxs_new) = if let Some(ref base) = ctx.base_xs {
transmission::broadened_cross_sections_with_derivative_from_base(
&ctx.energies,
base,
&ctx.resonance_data,
t_current,
ctx.instrument.as_ref(),
)
} else {
transmission::broadened_cross_sections_with_derivative(
&ctx.energies,
&ctx.resonance_data,
t_current,
ctx.instrument.as_ref(),
)
}
.map_err(|e| {
FittingError::EvaluationFailed(format!(
"poisson_fit_lbfgsb: broadening failed: {e}"
))
})?;
xs_owned = xs_new;
dxs_dt = dxs_new;
last_temperature = t_current;
}
}
let xs_ref: &[Vec<f64>] = if temp_ctx.is_some() {
&xs_owned
} else {
cross_sections
};
let y_model_now = model.evaluate(&all_vals_buf)?;
let t_now = compute_transmission(n_e, xs_ref, density_indices, &all_vals_buf);
let grad = compute_analytic_gradient(&AnalyticGradientCtx {
y_obs,
y_model: &y_model_now,
flux,
t_now: &t_now,
xs: xs_ref,
param_isotopes: ¶m_isotopes,
temp_free_pos,
density_indices,
all_vals: &all_vals_buf,
dxs_dt: &dxs_dt,
});
params.free_values_into(&mut free_vals_buf);
if let Some(ref pg) = prev_grad {
let s_vec: Vec<f64> = free_vals_buf
.iter()
.zip(prev_free.iter())
.map(|(&c, &p)| c - p)
.collect();
let y_vec: Vec<f64> = grad
.iter()
.zip(pg.iter())
.map(|(&gc, &gp)| gc - gp)
.collect();
lbfgs.push(s_vec, y_vec);
}
let zero_active_bounds = |vec: &mut [f64], g: &[f64]| {
for (j, &fi) in free_indices.iter().enumerate() {
let p = ¶ms.params[fi];
let at_lower = (p.value - p.lower).abs() < PIVOT_FLOOR && g[j] > 0.0;
let at_upper = (p.upper - p.value).abs() < PIVOT_FLOOR && g[j] < 0.0;
if at_lower || at_upper {
vec[j] = 0.0;
}
}
};
let mut projected_grad = grad.clone();
zero_active_bounds(&mut projected_grad, &grad);
let proj_grad_norm: f64 = projected_grad.iter().map(|g| g * g).sum::<f64>().sqrt();
if proj_grad_norm < config.tol_param {
converged = true;
break;
}
let mut direction = lbfgs.apply_two_loop(&projected_grad);
zero_active_bounds(&mut direction, &grad);
if !is_descent_direction(&direction, &grad) {
direction = projected_grad.clone();
}
let search_norm: f64 = direction.iter().map(|d| d * d).sum::<f64>().sqrt();
let initial_alpha = config.step_size / search_norm.max(1.0);
prev_free.clear();
prev_free.extend_from_slice(&free_vals_buf);
prev_grad = Some(grad.clone());
match backtracking_line_search(
model,
params,
y_obs,
&prev_free,
&direction,
initial_alpha,
config,
&grad,
nll,
&mut all_vals_buf,
&mut free_vals_buf,
&mut trial_free_buf,
) {
Some(new_nll) => nll = new_nll,
None => {
break;
}
}
params.free_values_into(&mut free_vals_buf);
let step_norm: f64 = prev_free
.iter()
.zip(free_vals_buf.iter())
.map(|(o, n)| (o - n).powi(2))
.sum::<f64>()
.sqrt();
if step_norm < config.tol_param {
converged = true;
break;
}
}
Ok(PoissonResult {
nll,
iterations: iter,
converged,
params: params.all_values(),
})
}
#[cfg(test)]
mod tests {
use super::*;
use crate::lm::FitModel;
use crate::parameters::FitParameter;
struct ExponentialModel {
x: Vec<f64>,
flux: Vec<f64>,
}
impl FitModel for ExponentialModel {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let b = params[0]; Ok(self
.x
.iter()
.zip(self.flux.iter())
.map(|(&xi, &fi)| fi * (-b * xi).exp())
.collect())
}
}
#[test]
fn test_poisson_nll_perfect_match() {
let y_obs = vec![10.0, 20.0, 30.0];
let y_model = vec![10.0, 20.0, 30.0];
let nll = poisson_nll(&y_obs, &y_model);
let expected: f64 = y_obs
.iter()
.zip(y_model.iter())
.map(|(&o, &m)| m - o * m.ln())
.sum();
assert!((nll - expected).abs() < 1e-10);
}
#[test]
fn test_poisson_fit_exponential() {
let x: Vec<f64> = (0..20).map(|i| i as f64 * 0.5).collect();
let true_b = 0.5;
let flux: Vec<f64> = vec![1000.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[true_b]).unwrap();
let mut params = ParameterSet::new(vec![
FitParameter::non_negative("b", 1.0), ]);
let result = poisson_fit(&model, &y_obs, &mut params, &PoissonConfig::default()).unwrap();
assert!(
result.converged,
"Poisson fit did not converge after {} iterations",
result.iterations,
);
assert!(
(result.params[0] - true_b).abs() / true_b < 0.05,
"Fitted b = {}, true = {}, error = {:.1}%",
result.params[0],
true_b,
(result.params[0] - true_b).abs() / true_b * 100.0,
);
}
#[test]
fn test_poisson_fit_low_counts() {
let x: Vec<f64> = (0..30).map(|i| i as f64 * 0.2).collect();
let true_b = 0.3;
let flux: Vec<f64> = vec![10.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[true_b]).unwrap();
let mut params = ParameterSet::new(vec![FitParameter::non_negative("b", 0.1)]);
let result = poisson_fit(&model, &y_obs, &mut params, &PoissonConfig::default()).unwrap();
assert!(result.converged);
assert!(
(result.params[0] - true_b).abs() / true_b < 0.1,
"Low-count: fitted b = {}, true = {}",
result.params[0],
true_b,
);
}
#[test]
fn test_poisson_non_negativity() {
let x: Vec<f64> = (0..10).map(|i| i as f64).collect();
let flux: Vec<f64> = vec![100.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs: Vec<f64> = vec![100.0; x.len()];
let mut params = ParameterSet::new(vec![FitParameter::non_negative("b", 1.0)]);
let result = poisson_fit(&model, &y_obs, &mut params, &PoissonConfig::default()).unwrap();
assert!(
result.params[0] >= 0.0,
"b should be non-negative, got {}",
result.params[0],
);
assert!(
result.params[0] < 0.1,
"b should be ~0, got {}",
result.params[0],
);
}
#[test]
fn test_counts_model() {
struct ConstTransmission;
impl FitModel for ConstTransmission {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(vec![params[0]; 3])
}
}
let t_model = ConstTransmission;
let flux = [100.0, 200.0, 300.0];
let background = [5.0, 10.0, 15.0];
let counts_model = CountsModel {
transmission_model: &t_model,
flux: &flux,
background: &background,
};
let result = counts_model.evaluate(&[0.5]).unwrap();
assert!((result[0] - 55.0).abs() < 1e-10);
assert!((result[1] - 110.0).abs() < 1e-10);
assert!((result[2] - 165.0).abs() < 1e-10);
}
#[test]
fn test_analytic_matches_fd_result() {
let x: Vec<f64> = (0..20).map(|i| i as f64 * 0.5).collect();
let true_b = 0.5;
let flux: Vec<f64> = vec![1000.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[true_b]).unwrap();
let cross_sections = vec![x.clone()];
let mut params_fd = ParameterSet::new(vec![FitParameter::non_negative("b", 1.0)]);
let res_fd =
poisson_fit(&model, &y_obs, &mut params_fd, &PoissonConfig::default()).unwrap();
let mut params_an = ParameterSet::new(vec![FitParameter::non_negative("b", 1.0)]);
let res_an = poisson_fit_analytic(
&model,
&y_obs,
&flux,
&cross_sections,
&[0],
&mut params_an,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(res_fd.converged, "FD did not converge");
assert!(res_an.converged, "Analytic did not converge");
assert!(
(res_fd.params[0] - res_an.params[0]).abs() < 1e-4,
"FD={}, Analytic={} should agree",
res_fd.params[0],
res_an.params[0],
);
assert!(
res_an.iterations <= res_fd.iterations,
"Analytic ({} iters) should not need more iterations than FD ({})",
res_an.iterations,
res_fd.iterations,
);
}
#[test]
fn test_analytic_two_isotopes() {
let n_e = 30;
let sigma1: Vec<f64> = (0..n_e).map(|i| 1.0 + 0.1 * i as f64).collect();
let sigma2: Vec<f64> = (0..n_e).map(|i| 0.5 + 0.05 * (n_e - i) as f64).collect();
let flux: Vec<f64> = vec![500.0; n_e];
let true_n1 = 0.3;
let true_n2 = 0.7;
struct TwoIsotopeModel {
sigma1: Vec<f64>,
sigma2: Vec<f64>,
flux: Vec<f64>,
}
impl FitModel for TwoIsotopeModel {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let (n1, n2) = (params[0], params[1]);
Ok(self
.sigma1
.iter()
.zip(self.sigma2.iter())
.zip(self.flux.iter())
.map(|((&s1, &s2), &f)| f * (-n1 * s1 - n2 * s2).exp())
.collect())
}
}
let model = TwoIsotopeModel {
sigma1: sigma1.clone(),
sigma2: sigma2.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[true_n1, true_n2]).unwrap();
let mut params = ParameterSet::new(vec![
FitParameter::non_negative("n1", 0.1),
FitParameter::non_negative("n2", 0.1),
]);
let result = poisson_fit_analytic(
&model,
&y_obs,
&flux,
&[sigma1, sigma2],
&[0, 1],
&mut params,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(result.converged, "Two-isotope analytic did not converge");
assert!(
(result.params[0] - true_n1).abs() / true_n1 < 0.01,
"n1: fitted={}, true={}",
result.params[0],
true_n1,
);
assert!(
(result.params[1] - true_n2).abs() / true_n2 < 0.01,
"n2: fitted={}, true={}",
result.params[1],
true_n2,
);
}
#[test]
fn test_analytic_zero_density_convergence() {
let x: Vec<f64> = (0..20).map(|i| i as f64 * 0.5).collect();
let flux: Vec<f64> = vec![1000.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[0.0]).unwrap();
let cross_sections = vec![x.clone()];
let mut params = ParameterSet::new(vec![FitParameter::non_negative("b", 0.5)]);
let result = poisson_fit_analytic(
&model,
&y_obs,
&flux,
&cross_sections,
&[0],
&mut params,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(result.converged, "Zero-density fit did not converge");
assert!(
result.params[0] < 0.01,
"b should be ~0, got {}",
result.params[0],
);
assert!(
result.params[0] >= 0.0,
"b must be non-negative, got {}",
result.params[0],
);
}
#[test]
fn test_analytic_nonzero_background() {
let n_e = 20;
let sigma: Vec<f64> = (0..n_e).map(|i| 1.0 + 0.2 * i as f64).collect();
let flux: Vec<f64> = vec![500.0; n_e];
let background: Vec<f64> = vec![20.0; n_e]; let true_b = 0.4;
struct PureTransmission {
sigma: Vec<f64>,
}
impl FitModel for PureTransmission {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let b = params[0];
Ok(self.sigma.iter().map(|&s| (-b * s).exp()).collect())
}
}
let t_model = PureTransmission {
sigma: sigma.clone(),
};
let counts_model = CountsModel {
transmission_model: &t_model,
flux: &flux,
background: &background,
};
let y_obs = counts_model.evaluate(&[true_b]).unwrap();
let mut params = ParameterSet::new(vec![FitParameter::non_negative("b", 0.2)]);
let config = PoissonConfig {
max_iter: 500,
..PoissonConfig::default()
};
let result = poisson_fit_analytic(
&counts_model,
&y_obs,
&flux,
&[sigma],
&[0],
&mut params,
&config,
None,
)
.unwrap();
assert!(
result.converged,
"Nonzero-background fit did not converge after {} iters",
result.iterations,
);
assert!(
(result.params[0] - true_b).abs() / true_b < 0.01,
"With B≠0: fitted={}, true={}, error={:.2}%",
result.params[0],
true_b,
(result.params[0] - true_b).abs() / true_b * 100.0,
);
}
#[test]
fn test_poisson_analytic_temperature_recovery() {
use nereids_core::types::Isotope;
use nereids_endf::resonance::{LGroup, Resonance, ResonanceFormalism, ResonanceRange};
use nereids_physics::transmission;
let resonance_data = vec![nereids_endf::resonance::ResonanceData {
isotope: Isotope::new(92, 238).unwrap(),
za: 92238,
awr: 236.006,
ranges: vec![ResonanceRange {
energy_low: 1e-5,
energy_high: 1e4,
resolved: true,
formalism: ResonanceFormalism::ReichMoore,
target_spin: 0.0,
scattering_radius: 9.4285,
naps: 1,
l_groups: vec![LGroup {
l: 0,
awr: 236.006,
apl: 0.0,
qx: 0.0,
lrx: 0,
resonances: vec![Resonance {
energy: 6.674,
j: 0.5,
gn: 1.493e-3,
gg: 23.0e-3,
gfa: 0.0,
gfb: 0.0,
}],
}],
rml: None,
urr: None,
ap_table: None,
r_external: vec![],
}],
}];
let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
let n_e = energies.len();
let true_density = 0.0005;
let true_temp = 300.0;
let flux: Vec<f64> = vec![10000.0; n_e];
let xs_true = transmission::broadened_cross_sections(
&energies,
&resonance_data,
true_temp,
None,
None,
)
.unwrap();
let y_obs: Vec<f64> = (0..n_e)
.map(|e| {
let t = (-true_density * xs_true[0][e]).exp();
flux[e] * t
})
.collect();
use nereids_physics::transmission::SampleParams;
struct PhysicsCountsModel {
energies: Vec<f64>,
resonance_data: Vec<nereids_endf::resonance::ResonanceData>,
flux: Vec<f64>,
}
impl FitModel for PhysicsCountsModel {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let density = params[0];
let temperature = params[1];
let sample =
SampleParams::new(temperature, vec![(self.resonance_data[0].clone(), density)])
.map_err(|e| FittingError::EvaluationFailed(e.to_string()))?;
let transmission = transmission::forward_model(&self.energies, &sample, None)
.map_err(|e| FittingError::EvaluationFailed(e.to_string()))?;
Ok(transmission
.iter()
.zip(self.flux.iter())
.map(|(&t, &f)| f * t)
.collect())
}
}
let model = PhysicsCountsModel {
energies: energies.clone(),
resonance_data: resonance_data.clone(),
flux: flux.clone(),
};
let init_temp = 200.0;
let xs_init = transmission::broadened_cross_sections(
&energies,
&resonance_data,
init_temp,
None,
None,
)
.unwrap();
let mut params = ParameterSet::new(vec![
FitParameter::non_negative("density", 0.001),
FitParameter {
name: "temperature".into(),
value: init_temp,
lower: 50.0,
upper: 1000.0,
fixed: false,
},
]);
let temp_ctx = TemperatureContext {
temperature_index: 1,
resonance_data: resonance_data.clone(),
energies: energies.clone(),
instrument: None,
base_xs: None,
};
let config = PoissonConfig {
max_iter: 5000,
tol_param: 1e-12,
..PoissonConfig::default()
};
let result = poisson_fit_analytic(
&model,
&y_obs,
&flux,
&xs_init,
&[0],
&mut params,
&config,
Some(&temp_ctx),
)
.unwrap();
assert!(
result.converged,
"Temperature+density fit did not converge after {} iterations",
result.iterations,
);
let fitted_density = result.params[0];
let fitted_temp = result.params[1];
assert!(
(fitted_density - true_density).abs() / true_density < 0.01,
"density: fitted={}, true={}, error={:.2}%",
fitted_density,
true_density,
(fitted_density - true_density).abs() / true_density * 100.0,
);
assert!(
(fitted_temp - true_temp).abs() / true_temp < 0.01,
"temperature: fitted={:.1}, true={:.1}, error={:.2}%",
fitted_temp,
true_temp,
(fitted_temp - true_temp).abs() / true_temp * 100.0,
);
}
#[test]
fn test_poisson_nll_term_negative_model() {
let nll = poisson_nll_term(10.0, -5.0);
assert!(
nll.is_finite(),
"NLL should be finite for negative model, got {nll}"
);
let nll_at_eps = poisson_nll_term(10.0, POISSON_EPSILON);
assert!(
nll > nll_at_eps,
"NLL at mdl=-5 ({nll}) should exceed NLL at epsilon ({nll_at_eps})"
);
}
#[test]
fn test_poisson_nll_term_zero_obs() {
let nll = poisson_nll_term(0.0, 10.0);
assert!(
(nll - 10.0).abs() < 1e-10,
"NLL(obs=0, mdl=10) should be 10.0, got {nll}"
);
let nll_neg = poisson_nll_term(0.0, -1.0);
let nll_zero = poisson_nll_term(0.0, 0.0);
assert!(
nll_neg > nll_zero,
"NLL should grow as mdl goes more negative: nll(-1)={nll_neg} vs nll(0)={nll_zero}"
);
}
#[test]
fn test_all_fixed_params_nan_model_poisson() {
struct NanModel;
impl FitModel for NanModel {
fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(vec![f64::NAN; 5])
}
}
let y_obs = vec![10.0; 5];
let mut params = ParameterSet::new(vec![FitParameter::fixed("a", 1.0)]);
let result =
poisson_fit(&NanModel, &y_obs, &mut params, &PoissonConfig::default()).unwrap();
assert!(
!result.converged,
"All-fixed NaN model should not converge in Poisson fit"
);
assert!(!result.nll.is_finite(), "NLL should be non-finite");
assert_eq!(result.iterations, 0);
}
#[test]
fn test_all_fixed_params_nan_model_poisson_analytic() {
struct NanModel;
impl FitModel for NanModel {
fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(vec![f64::NAN; 5])
}
}
let y_obs = vec![10.0; 5];
let flux = vec![1000.0; 5];
let cross_sections = vec![vec![1.0; 5]];
let density_indices = vec![0];
let mut params = ParameterSet::new(vec![FitParameter::fixed("density", 0.5)]);
let result = poisson_fit_analytic(
&NanModel,
&y_obs,
&flux,
&cross_sections,
&density_indices,
&mut params,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(
!result.converged,
"All-fixed NaN model should not converge in poisson_fit_analytic"
);
assert!(!result.nll.is_finite(), "NLL should be non-finite");
assert_eq!(result.iterations, 0);
}
#[test]
fn test_all_fixed_params_good_model_poisson() {
let x: Vec<f64> = (0..10).map(|i| i as f64 * 0.5).collect();
let flux: Vec<f64> = vec![100.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[0.5]).unwrap();
let mut params = ParameterSet::new(vec![FitParameter::fixed("b", 0.5)]);
let result = poisson_fit(&model, &y_obs, &mut params, &PoissonConfig::default()).unwrap();
assert!(
result.converged,
"All-fixed good model should converge in Poisson fit"
);
assert!(result.nll.is_finite(), "NLL should be finite");
assert_eq!(result.iterations, 0);
}
#[test]
fn test_all_fixed_params_nan_model_poisson_lbfgsb() {
struct NanModel;
impl FitModel for NanModel {
fn evaluate(&self, _params: &[f64]) -> Result<Vec<f64>, FittingError> {
Ok(vec![f64::NAN; 5])
}
}
let y_obs = vec![10.0; 5];
let flux = vec![1000.0; 5];
let cross_sections = vec![vec![1.0; 5]];
let density_indices = vec![0];
let mut params = ParameterSet::new(vec![FitParameter::fixed("density", 0.5)]);
let result = poisson_fit_lbfgsb(
&NanModel,
&y_obs,
&flux,
&cross_sections,
&density_indices,
&mut params,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(
!result.converged,
"All-fixed NaN model should not converge in poisson_fit_lbfgsb"
);
assert!(!result.nll.is_finite(), "NLL should be non-finite");
assert_eq!(result.iterations, 0);
}
#[test]
fn test_lbfgsb_matches_analytic_result() {
let x: Vec<f64> = (0..20).map(|i| i as f64 * 0.5).collect();
let true_b = 0.5;
let flux: Vec<f64> = vec![1000.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[true_b]).unwrap();
let cross_sections = vec![x.clone()];
let mut params_an = ParameterSet::new(vec![FitParameter::non_negative("b", 1.0)]);
let res_an = poisson_fit_analytic(
&model,
&y_obs,
&flux,
&cross_sections,
&[0],
&mut params_an,
&PoissonConfig::default(),
None,
)
.unwrap();
let mut params_lb = ParameterSet::new(vec![FitParameter::non_negative("b", 1.0)]);
let res_lb = poisson_fit_lbfgsb(
&model,
&y_obs,
&flux,
&cross_sections,
&[0],
&mut params_lb,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(res_an.converged, "Analytic did not converge");
assert!(res_lb.converged, "L-BFGS-B did not converge");
assert!(
(res_an.params[0] - res_lb.params[0]).abs() < 1e-4,
"Analytic={}, L-BFGS-B={} should agree",
res_an.params[0],
res_lb.params[0],
);
}
#[test]
fn test_lbfgsb_two_isotopes() {
let n_e = 30;
let sigma1: Vec<f64> = (0..n_e).map(|i| 1.0 + 0.1 * i as f64).collect();
let sigma2: Vec<f64> = (0..n_e).map(|i| 0.5 + 0.05 * (n_e - i) as f64).collect();
let flux: Vec<f64> = vec![500.0; n_e];
let true_n1 = 0.3;
let true_n2 = 0.7;
struct TwoIsotopeModel {
sigma1: Vec<f64>,
sigma2: Vec<f64>,
flux: Vec<f64>,
}
impl FitModel for TwoIsotopeModel {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let (n1, n2) = (params[0], params[1]);
Ok(self
.sigma1
.iter()
.zip(self.sigma2.iter())
.zip(self.flux.iter())
.map(|((&s1, &s2), &f)| f * (-n1 * s1 - n2 * s2).exp())
.collect())
}
}
let model = TwoIsotopeModel {
sigma1: sigma1.clone(),
sigma2: sigma2.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[true_n1, true_n2]).unwrap();
let mut params = ParameterSet::new(vec![
FitParameter::non_negative("n1", 0.1),
FitParameter::non_negative("n2", 0.1),
]);
let result = poisson_fit_lbfgsb(
&model,
&y_obs,
&flux,
&[sigma1, sigma2],
&[0, 1],
&mut params,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(
result.converged,
"L-BFGS-B two-isotope did not converge after {} iters",
result.iterations,
);
assert!(
(result.params[0] - true_n1).abs() / true_n1 < 0.01,
"n1: fitted={}, true={}",
result.params[0],
true_n1,
);
assert!(
(result.params[1] - true_n2).abs() / true_n2 < 0.01,
"n2: fitted={}, true={}",
result.params[1],
true_n2,
);
}
#[test]
fn test_lbfgsb_temperature_recovery() {
use nereids_core::types::Isotope;
use nereids_endf::resonance::{LGroup, Resonance, ResonanceFormalism, ResonanceRange};
use nereids_physics::transmission;
let resonance_data = vec![nereids_endf::resonance::ResonanceData {
isotope: Isotope::new(92, 238).unwrap(),
za: 92238,
awr: 236.006,
ranges: vec![ResonanceRange {
energy_low: 1e-5,
energy_high: 1e4,
resolved: true,
formalism: ResonanceFormalism::ReichMoore,
target_spin: 0.0,
scattering_radius: 9.4285,
naps: 1,
l_groups: vec![LGroup {
l: 0,
awr: 236.006,
apl: 0.0,
qx: 0.0,
lrx: 0,
resonances: vec![Resonance {
energy: 6.674,
j: 0.5,
gn: 1.493e-3,
gg: 23.0e-3,
gfa: 0.0,
gfb: 0.0,
}],
}],
rml: None,
urr: None,
ap_table: None,
r_external: vec![],
}],
}];
let energies: Vec<f64> = (0..201).map(|i| 4.0 + (i as f64) * 0.025).collect();
let n_e = energies.len();
let true_density = 0.0005;
let true_temp = 300.0;
let flux: Vec<f64> = vec![10000.0; n_e];
let xs_true = transmission::broadened_cross_sections(
&energies,
&resonance_data,
true_temp,
None,
None,
)
.unwrap();
let y_obs: Vec<f64> = (0..n_e)
.map(|e| {
let t = (-true_density * xs_true[0][e]).exp();
flux[e] * t
})
.collect();
use nereids_physics::transmission::SampleParams;
struct PhysicsCountsModel {
energies: Vec<f64>,
resonance_data: Vec<nereids_endf::resonance::ResonanceData>,
flux: Vec<f64>,
}
impl FitModel for PhysicsCountsModel {
fn evaluate(&self, params: &[f64]) -> Result<Vec<f64>, FittingError> {
let density = params[0];
let temperature = params[1];
let sample =
SampleParams::new(temperature, vec![(self.resonance_data[0].clone(), density)])
.map_err(|e| FittingError::EvaluationFailed(e.to_string()))?;
let transmission = transmission::forward_model(&self.energies, &sample, None)
.map_err(|e| FittingError::EvaluationFailed(e.to_string()))?;
Ok(transmission
.iter()
.zip(self.flux.iter())
.map(|(&t, &f)| f * t)
.collect())
}
}
let model = PhysicsCountsModel {
energies: energies.clone(),
resonance_data: resonance_data.clone(),
flux: flux.clone(),
};
let init_temp = 200.0;
let xs_init = transmission::broadened_cross_sections(
&energies,
&resonance_data,
init_temp,
None,
None,
)
.unwrap();
let mut params = ParameterSet::new(vec![
FitParameter::non_negative("density", 0.001),
FitParameter {
name: "temperature".into(),
value: init_temp,
lower: 50.0,
upper: 1000.0,
fixed: false,
},
]);
let temp_ctx = TemperatureContext {
temperature_index: 1,
resonance_data: resonance_data.clone(),
energies: energies.clone(),
instrument: None,
base_xs: None,
};
let config = PoissonConfig {
max_iter: 200,
tol_param: 1e-12,
..PoissonConfig::default()
};
let result = poisson_fit_lbfgsb(
&model,
&y_obs,
&flux,
&xs_init,
&[0],
&mut params,
&config,
Some(&temp_ctx),
)
.unwrap();
assert!(
result.converged,
"L-BFGS-B temperature+density fit did not converge after {} iterations",
result.iterations,
);
let fitted_density = result.params[0];
let fitted_temp = result.params[1];
assert!(
(fitted_density - true_density).abs() / true_density < 0.01,
"density: fitted={}, true={}, error={:.2}%",
fitted_density,
true_density,
(fitted_density - true_density).abs() / true_density * 100.0,
);
assert!(
(fitted_temp - true_temp).abs() / true_temp < 0.01,
"temperature: fitted={:.1}, true={:.1}, error={:.2}%",
fitted_temp,
true_temp,
(fitted_temp - true_temp).abs() / true_temp * 100.0,
);
}
#[test]
fn test_lbfgsb_zero_density() {
let x: Vec<f64> = (0..20).map(|i| i as f64 * 0.5).collect();
let flux: Vec<f64> = vec![1000.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[0.0]).unwrap();
let cross_sections = vec![x.clone()];
let mut params = ParameterSet::new(vec![FitParameter::non_negative("b", 0.5)]);
let result = poisson_fit_lbfgsb(
&model,
&y_obs,
&flux,
&cross_sections,
&[0],
&mut params,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(
result.converged,
"L-BFGS-B zero-density fit did not converge"
);
assert!(
result.params[0] < 0.01,
"b should be ~0, got {}",
result.params[0],
);
assert!(
result.params[0] >= 0.0,
"b must be non-negative, got {}",
result.params[0],
);
}
#[test]
fn test_lbfgsb_all_fixed() {
let x: Vec<f64> = (0..10).map(|i| i as f64 * 0.5).collect();
let flux: Vec<f64> = vec![100.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[0.5]).unwrap();
let cross_sections = vec![x.clone()];
let mut params = ParameterSet::new(vec![FitParameter::fixed("b", 0.5)]);
let result = poisson_fit_lbfgsb(
&model,
&y_obs,
&flux,
&cross_sections,
&[0],
&mut params,
&PoissonConfig::default(),
None,
)
.unwrap();
assert!(
result.converged,
"All-fixed L-BFGS-B should converge immediately"
);
assert!(result.nll.is_finite(), "NLL should be finite");
assert_eq!(result.iterations, 0);
}
#[test]
fn test_lbfgsb_memory_one() {
let x: Vec<f64> = (0..20).map(|i| i as f64 * 0.5).collect();
let true_b = 0.5;
let flux: Vec<f64> = vec![1000.0; x.len()];
let model = ExponentialModel {
x: x.clone(),
flux: flux.clone(),
};
let y_obs = model.evaluate(&[true_b]).unwrap();
let cross_sections = vec![x.clone()];
let mut params = ParameterSet::new(vec![FitParameter::non_negative("b", 1.0)]);
let config = PoissonConfig {
lbfgsb_memory: 1,
max_iter: 500, ..PoissonConfig::default()
};
let result = poisson_fit_lbfgsb(
&model,
&y_obs,
&flux,
&cross_sections,
&[0],
&mut params,
&config,
None,
)
.unwrap();
assert!(
result.converged,
"L-BFGS-B with m=1 did not converge after {} iterations",
result.iterations,
);
assert!(
(result.params[0] - true_b).abs() / true_b < 0.01,
"m=1: fitted={}, true={}, error={:.2}%",
result.params[0],
true_b,
(result.params[0] - true_b).abs() / true_b * 100.0,
);
}
#[test]
fn test_is_descent_direction_opposite_vectors() {
assert!(is_descent_direction(&[1.0, 0.0], &[1.0, 0.0]));
assert!(is_descent_direction(&[0.5, 0.3], &[2.0, 1.0]));
}
#[test]
fn test_is_descent_direction_orthogonal_not_descent() {
assert!(!is_descent_direction(&[0.0, 1.0], &[1.0, 0.0]));
assert!(!is_descent_direction(&[1.0, 0.0], &[0.0, 1.0]));
}
#[test]
fn test_is_descent_direction_ascending_not_descent() {
assert!(!is_descent_direction(&[-1.0, 0.0], &[1.0, 0.0]));
assert!(!is_descent_direction(&[-0.5, -0.3], &[2.0, 1.0]));
}
}