lcdm_core/params.rs
1// lcdm-core/src/params.rs
2
3use crate::constants::{C, G, H0_UNITS_TO_INV_SEC, SIGMA_SB};
4#[cfg(feature = "serde")]
5use serde::{Deserialize, Serialize};
6
7use crate::spec::{AccuracyPreset, DerivedParams, NeutrinoSpec};
8
9#[derive(Debug, Clone)]
10#[cfg_attr(feature = "serde", derive(Serialize, Deserialize))]
11pub struct ScientificParams {
12 /// Canonical, user-facing parameter set (inputs in convenient forms).
13 pub canonical: CanonicalParams,
14
15 /// Explicit neutrino model specification (effective vs split).
16 pub spec: NeutrinoSpec,
17
18 /// Derived quantities computed from the canonical/spec inputs.
19 pub derived: DerivedParams,
20
21 /// Global accuracy / performance preset controlling numerical knobs.
22 pub accuracy: AccuracyPreset,
23}
24
25#[derive(Debug, Clone)]
26#[cfg_attr(feature = "serde", derive(Serialize, Deserialize))]
27pub struct CanonicalParams {
28 /// Dimensionless Hubble parameter h, where H0 = 100 h (km/s)/Mpc.
29 pub h: f64,
30
31 /// CMB temperature today in Kelvin.
32 pub T_cmb0: f64,
33
34 /// Physical baryon density: ω_b ≡ Ω_b h².
35 pub omega_b: f64,
36
37 /// Physical cold dark matter density: ω_cdm ≡ Ω_cdm h².
38 pub omega_cdm: f64,
39
40 /// Physical dark energy density: ω_Λ ≡ Ω_Λ h².
41 pub omega_lambda: f64,
42
43 /// Physical curvature density: ω_k ≡ Ω_k h².
44 pub omega_k: f64,
45
46 /// Effective number of relativistic species (or equivalent accounting in split mode).
47 pub n_eff: f64,
48
49 /// Neutrino masses in eV (interpretation depends on the neutrino spec).
50 pub m_nu: Vec<f64>,
51
52 /// Helium mass fraction (placeholder / future BBN integration).
53 pub Y_He: f64,
54
55 /// Dark energy equation-of-state parameter w0 (CPL form).
56 pub w0: f64,
57
58 /// Dark energy equation-of-state running parameter wa (CPL form).
59 pub wa: f64,
60}
61
62impl CanonicalParams {
63 /// Hubble constant H0 in (km/s)/Mpc.
64 pub fn H0(&self) -> f64 {
65 self.h * 100.0
66 }
67}
68
69#[derive(Debug, Clone, Default)]
70pub struct CosmoBuilder {
71 h: Option<f64>,
72 h0: Option<f64>,
73
74 // Physical densities (ω = Ω h^2).
75 omega_b_phys: Option<f64>,
76 omega_cdm_phys: Option<f64>,
77
78 // Dimensionless densities (Ω).
79 Omega_b: Option<f64>,
80 Omega_cdm: Option<f64>,
81 Omega_m: Option<f64>,
82 Omega_k: Option<f64>,
83 Omega_L: Option<f64>,
84
85 T_cmb0: Option<f64>,
86
87 // Explicit neutrino model specification.
88 nu_model: Option<NeutrinoSpec>,
89
90 w0: Option<f64>,
91 wa: Option<f64>,
92 accuracy: AccuracyPreset,
93 flat: bool,
94}
95
96impl CosmoBuilder {
97 /// Create a new builder with default accuracy settings.
98 pub fn new() -> Self {
99 Self {
100 accuracy: AccuracyPreset::Default,
101 ..Self::default()
102 }
103 }
104
105 /// Enforce spatial flatness (Ω_k = 0) and solve for ω_Λ by closure.
106 pub fn flat(mut self) -> Self {
107 self.flat = true;
108 self
109 }
110
111 /// Set the dimensionless Hubble parameter h (H0 = 100 h).
112 pub fn h(mut self, val: f64) -> Self {
113 self.h = Some(val);
114 self
115 }
116
117 /// Set H0 in (km/s)/Mpc (converted internally to h).
118 pub fn H0(mut self, val: f64) -> Self {
119 self.h0 = Some(val);
120 self
121 }
122
123 // -------------------------------------------------------------------------
124 // Dimensionless density API (Ω)
125 // -------------------------------------------------------------------------
126
127 #[allow(non_snake_case)]
128 /// Set total matter density Ω_m.
129 pub fn Omega_m(mut self, val: f64) -> Self {
130 self.Omega_m = Some(val);
131 self
132 }
133
134 #[allow(non_snake_case)]
135 /// Set baryon density Ω_b.
136 pub fn Omega_b(mut self, val: f64) -> Self {
137 self.Omega_b = Some(val);
138 self
139 }
140
141 #[allow(non_snake_case)]
142 /// Set CDM density Ω_cdm.
143 pub fn Omega_cdm(mut self, val: f64) -> Self {
144 self.Omega_cdm = Some(val);
145 self
146 }
147
148 #[allow(non_snake_case)]
149 /// Set curvature density Ω_k.
150 pub fn Omega_k(mut self, val: f64) -> Self {
151 self.Omega_k = Some(val);
152 self
153 }
154
155 #[allow(non_snake_case)]
156 /// Set dark energy density Ω_Λ.
157 pub fn Omega_L(mut self, val: f64) -> Self {
158 self.Omega_L = Some(val);
159 self
160 }
161
162 // -------------------------------------------------------------------------
163 // Physical density API (ω ≡ Ω h²)
164 // -------------------------------------------------------------------------
165
166 /// Set physical baryon density ω_b.
167 pub fn omega_b_phys(mut self, val: f64) -> Self {
168 self.omega_b_phys = Some(val);
169 self
170 }
171
172 /// Set physical CDM density ω_cdm.
173 pub fn omega_cdm_phys(mut self, val: f64) -> Self {
174 self.omega_cdm_phys = Some(val);
175 self
176 }
177
178 // -------------------------------------------------------------------------
179 // Misc configuration
180 // -------------------------------------------------------------------------
181
182 #[allow(non_snake_case)]
183 /// Set the present-day CMB temperature T_cmb0 in Kelvin.
184 pub fn T_cmb0(mut self, val: f64) -> Self {
185 self.T_cmb0 = Some(val);
186 self
187 }
188
189 /// Set the global accuracy preset controlling numerical settings.
190 pub fn accuracy(mut self, val: AccuracyPreset) -> Self {
191 self.accuracy = val;
192 self
193 }
194
195 // -------------------------------------------------------------------------
196 // Neutrino model selection
197 // -------------------------------------------------------------------------
198
199 /// Use an "effective" neutrino model with total `n_eff` and a mass list.
200 ///
201 /// For massless neutrinos, pass `[0.0]` (or multiple zeros if you want explicit states).
202 pub fn neutrino_effective(mut self, n_eff: f64, masses_ev: Vec<f64>) -> Self {
203 self.nu_model = Some(NeutrinoSpec::Effective { n_eff, masses_ev });
204 self
205 }
206
207 /// Use a "split" neutrino model with explicit ultra-relativistic contribution `n_ur`
208 /// plus a list of massive eigenstates and per-species temperature factors.
209 pub fn neutrino_split(mut self, n_ur: f64, masses_ev: Vec<f64>, temp_factors: Vec<f64>) -> Self {
210 self.nu_model = Some(NeutrinoSpec::Split {
211 n_ur,
212 masses_ev,
213 temp_factors,
214 });
215 self
216 }
217
218 /// Build `ScientificParams` by resolving inputs, enforcing consistency, and computing
219 /// all derived quantities.
220 pub fn build(self) -> Result<ScientificParams, String> {
221 // Resolve h from either h or H0 input (and validate if both were provided).
222 let h = match (self.h, self.h0) {
223 (Some(hv), None) => hv,
224 (None, Some(h0v)) => h0v / 100.0,
225 (Some(hv), Some(h0v)) => {
226 if (hv * 100.0 - h0v).abs() > 1e-5 {
227 return Err("Inconsistent H0/h".into());
228 }
229 hv
230 }
231 (None, None) => return Err("Specify h or H0".into()),
232 };
233
234 // Default to the commonly used CMB temperature if not specified.
235 let t_cmb0 = self.T_cmb0.unwrap_or(2.7255);
236
237 // Neutrino model is mandatory (explicit).
238 let spec = self.nu_model.ok_or("Neutrino model must be explicitly set")?;
239
240 // Photon energy density today: ρ_γ = (4σ/c) T^4 (from blackbody radiation).
241 let rho_g_energy = 4.0 * SIGMA_SB / C * t_cmb0.powi(4);
242
243 // Critical density for h=1, using H100 = 100 (km/s)/Mpc converted to SI.
244 let h100_si = 100.0 * H0_UNITS_TO_INV_SEC;
245 let rho_crit_100 = 3.0 * h100_si * h100_si / (8.0 * std::f64::consts::PI * G);
246
247 // Photon physical density ω_γ0 = ρ_γ / ρ_crit(h=1).
248 let omega_g = (rho_g_energy / (C * C)) / rho_crit_100;
249
250 // Base neutrino temperature: T_ν0 = T_γ0 (4/11)^(1/3).
251 let t_nu_base = t_cmb0 * (4.0 / 11.0_f64).powf(1.0 / 3.0);
252
253 // Numerical knobs for neutrino integration.
254 let y_max = self.accuracy.neutrino_y_max();
255 let steps = self.accuracy.neutrino_integ_steps();
256
257 // Resolve ω_ur0 (explicit UR component), ω_ν0 (integrated massive/effective),
258 // and the effective neutrino temperature used in the calculation.
259 let (omega_ur0, omega_nu0, t_nu_eff) = match &spec {
260 NeutrinoSpec::Effective { n_eff, masses_ev } => {
261 // For effective mode, masses must be provided (use [0.0] for massless).
262 if masses_ev.is_empty() {
263 return Err(
264 "Effective neutrino model: masses_ev must be non-empty (use [0.0] for massless)"
265 .into(),
266 );
267 }
268
269 // Distribute N_eff across the provided number of states via an effective temperature.
270 let n = masses_ev.len().max(1) as f64;
271 let t_eff = t_nu_base * (n_eff / n).powf(0.25);
272
273 // Hardcode degeneracy_g = 2.0 for each species.
274 let o_nu = crate::neutrino::calc_omega_nu_z_cfg(
275 masses_ev,
276 t_eff,
277 0.0,
278 y_max,
279 steps,
280 2.0,
281 );
282 (0.0, o_nu, t_eff)
283 }
284 NeutrinoSpec::Split {
285 n_ur,
286 masses_ev,
287 temp_factors,
288 } => {
289 // Validate split-model inputs.
290 if *n_ur < 0.0 {
291 return Err("n_ur < 0".into());
292 }
293 if masses_ev.len() != temp_factors.len() {
294 return Err("Length mismatch".into());
295 }
296 if masses_ev.iter().any(|&m| m < 1e-10) {
297 return Err("Split massive must be > 0".into());
298 }
299 if temp_factors.iter().any(|&tf| tf <= 0.0) {
300 return Err("Temp factors must be positive".into());
301 }
302
303 // Ultra-relativistic contribution uses the standard N_eff scaling:
304 // ρ_ν,UR / ρ_γ = (7/8) * (4/11)^(4/3) * N_eff
305 let o_ur =
306 *n_ur * (7.0 / 8.0) * (4.0 / 11.0_f64).powf(4.0 / 3.0) * omega_g;
307
308 // Integrate each massive eigenstate with its temperature factor.
309 let mut o_nu_total = 0.0;
310 for (&m, &tf) in masses_ev.iter().zip(temp_factors.iter()) {
311 o_nu_total += crate::neutrino::calc_omega_nu_z_cfg(
312 &[m],
313 t_nu_base * tf,
314 0.0,
315 y_max,
316 steps,
317 2.0,
318 );
319 }
320
321 (o_ur, o_nu_total, t_nu_base)
322 }
323 };
324
325 // ---------------------------------------------------------------------
326 // Resolve baryon density ω_b
327 // ---------------------------------------------------------------------
328
329 let omega_b = if let (Some(wb), Some(Ob)) = (self.omega_b_phys, self.Omega_b) {
330 // If both are provided, enforce consistency.
331 let wb2 = Ob * h * h;
332 if (wb - wb2).abs() > 1e-10 {
333 return Err("Inconsistent omega_b_phys vs Omega_b".into());
334 }
335 wb
336 } else if let Some(wb) = self.omega_b_phys {
337 wb
338 } else if let Some(Ob) = self.Omega_b {
339 Ob * h * h
340 } else {
341 return Err("Specify Omega_b or omega_b_phys".into());
342 };
343
344 // ---------------------------------------------------------------------
345 // Resolve CDM density ω_cdm
346 // ---------------------------------------------------------------------
347
348 let omega_cdm = if let (Some(wc), Some(Oc)) = (self.omega_cdm_phys, self.Omega_cdm) {
349 // If both are provided, enforce consistency.
350 let wc2 = Oc * h * h;
351 if (wc - wc2).abs() > 1e-10 {
352 return Err("Inconsistent omega_cdm_phys vs Omega_cdm".into());
353 }
354 wc
355 } else if let Some(wc) = self.omega_cdm_phys {
356 wc
357 } else if let Some(Oc) = self.Omega_cdm {
358 Oc * h * h
359 } else if let Some(Om) = self.Omega_m {
360 // If only Ω_m is given, infer ω_cdm = Ω_m h^2 - ω_b.
361 Om * h * h - omega_b
362 } else {
363 return Err("Specify Omega_cdm, omega_cdm_phys, or Omega_m".into());
364 };
365
366 // Prevent nonsensical negative CDM (e.g., Ω_m < Ω_b).
367 if omega_cdm < 0.0 {
368 return Err("Negative CDM (Omega_m < Omega_b)".into());
369 }
370
371 let omega_m = omega_b + omega_cdm;
372 let omega_non_de = omega_g + omega_ur0 + omega_nu0 + omega_m;
373
374 // ---------------------------------------------------------------------
375 // Resolve dark energy and curvature via closure
376 // ---------------------------------------------------------------------
377
378 let (omega_lambda, omega_k) = if self.flat {
379 // In flat mode, Ω_k must be 0 and ω_Λ is solved by closure:
380 // h^2 = ω_non_de + ω_Λ
381 if let Some(ok) = self.Omega_k {
382 if ok.abs() > 1e-12 {
383 return Err("flat() is set but Omega_k was provided (must be 0)".into());
384 }
385 }
386
387 let ol = h * h - omega_non_de;
388 if ol < 0.0 {
389 return Err("Total density > h^2 in flat model".into());
390 }
391
392 (ol, 0.0)
393 } else {
394 // Non-flat: allow explicit Ω_Λ and/or Ω_k if provided, otherwise solve by closure.
395 let ol = if let Some(o_l_in) = self.Omega_L {
396 let ol_sq = o_l_in * h * h;
397
398 // If Ω_k is also provided, ensure the residual is consistent.
399 if let Some(o_k_in) = self.Omega_k {
400 let ok_sq = o_k_in * h * h;
401 let residual = (h * h - omega_non_de - ol_sq - ok_sq).abs();
402 if residual > 1e-6 {
403 return Err("Inconsistent Omega_L and Omega_k".into());
404 }
405 }
406
407 ol_sq
408 } else {
409 // If Ω_Λ is not provided, infer it after subtracting any provided Ω_k.
410 h * h - omega_non_de - self.Omega_k.unwrap_or(0.0) * h * h
411 };
412
413 // ω_k is what remains after enforcing closure: h^2 = ω_non_de + ω_Λ + ω_k.
414 (ol, h * h - omega_non_de - ol)
415 };
416
417 // Build the derived-parameter cache.
418 let derived = DerivedParams::new(
419 h,
420 t_cmb0,
421 omega_b,
422 omega_cdm,
423 omega_k,
424 omega_lambda,
425 omega_g,
426 omega_ur0,
427 omega_nu0,
428 t_nu_base,
429 t_nu_eff,
430 );
431
432 // Canonical n_eff and m_nu bookkeeping:
433 // - Effective: pass through n_eff and the provided mass list.
434 // - Split: treat n_eff as (n_ur + number_of_massive_states).
435 let (n_eff_can, m_nu_can) = match &spec {
436 NeutrinoSpec::Effective { n_eff, masses_ev, .. } => (*n_eff, masses_ev.clone()),
437 NeutrinoSpec::Split { n_ur, masses_ev, .. } => (*n_ur + masses_ev.len() as f64, masses_ev.clone()),
438 };
439
440 Ok(ScientificParams {
441 canonical: CanonicalParams {
442 h,
443 T_cmb0: t_cmb0,
444 omega_b,
445 omega_cdm,
446 omega_lambda,
447 omega_k,
448 n_eff: n_eff_can,
449 m_nu: m_nu_can,
450 Y_He: 0.24,
451 w0: self.w0.unwrap_or(-1.0),
452 wa: self.wa.unwrap_or(0.0),
453 },
454 spec,
455 derived,
456 accuracy: self.accuracy,
457 })
458 }
459}