fusion_core/
particles.rs

1// SPDX-License-Identifier: AGPL-3.0-or-later | Commercial license available
2// © Concepts 1996–2026 Miroslav Šotek. All rights reserved.
3// © Code 2020–2026 Miroslav Šotek. All rights reserved.
4// ORCID: 0009-0009-3560-0851
5// Contact: www.anulum.li | protoscience@anulum.li
6//! Reduced particle tracker overlay for hybrid MHD/PIC-style current feedback.
7//!
8//! This module introduces a deterministic charged-particle pusher using the
9//! Boris integrator and a simple toroidal current deposition path on the
10//! Grad-Shafranov grid.
11
12use fusion_types::error::{FusionError, FusionResult};
13use fusion_types::state::Grid2D;
14use ndarray::{Array1, Array2};
15use std::f64::consts::PI;
16
17const MIN_RADIUS_M: f64 = 1e-9;
18const MIN_CELL_AREA_M2: f64 = 1e-12;
19const MIN_CURRENT_INTEGRAL: f64 = 1e-9;
20const ELEMENTARY_CHARGE_C: f64 = 1.602_176_634e-19;
21const ALPHA_MASS_KG: f64 = 6.644_657_335_7e-27;
22const ALPHA_CHARGE_C: f64 = 2.0 * ELEMENTARY_CHARGE_C;
23
24// Coulomb collision constants
25const ELECTRON_MASS_KG: f64 = 9.109_383_701_5e-31;
26const PROTON_MASS_KG: f64 = 1.672_621_923_69e-27;
27const VACUUM_PERMITTIVITY: f64 = 8.854_187_812_8e-12;
28const BOLTZMANN_J_PER_KEV: f64 = 1.602_176_634e-16;
29
30/// Charged macro-particle state.
31#[derive(Debug, Clone, Copy, PartialEq)]
32pub struct ChargedParticle {
33    pub x_m: f64,
34    pub y_m: f64,
35    pub z_m: f64,
36    pub vx_m_s: f64,
37    pub vy_m_s: f64,
38    pub vz_m_s: f64,
39    pub charge_c: f64,
40    pub mass_kg: f64,
41    pub weight: f64,
42}
43
44impl ChargedParticle {
45    /// Cylindrical major radius R from Cartesian position.
46    pub fn cylindrical_radius_m(&self) -> f64 {
47        (self.x_m * self.x_m + self.y_m * self.y_m).sqrt()
48    }
49
50    /// Toroidal velocity component v_phi in cylindrical basis.
51    pub fn toroidal_velocity_m_s(&self) -> f64 {
52        let r = self.cylindrical_radius_m().max(MIN_RADIUS_M);
53        (-self.y_m * self.vx_m_s + self.x_m * self.vy_m_s) / r
54    }
55
56    /// Non-relativistic kinetic energy [J].
57    pub fn kinetic_energy_j(&self) -> f64 {
58        let v2 = self.vx_m_s * self.vx_m_s + self.vy_m_s * self.vy_m_s + self.vz_m_s * self.vz_m_s;
59        0.5 * self.mass_kg * v2
60    }
61
62    /// Non-relativistic kinetic energy [MeV].
63    pub fn kinetic_energy_mev(&self) -> f64 {
64        self.kinetic_energy_j() / (1.0e6 * ELEMENTARY_CHARGE_C)
65    }
66}
67
68#[derive(Debug, Clone, Copy, PartialEq)]
69pub struct ParticlePopulationSummary {
70    pub count: usize,
71    pub mean_energy_mev: f64,
72    pub p95_energy_mev: f64,
73    pub max_energy_mev: f64,
74    pub runaway_fraction: f64,
75}
76
77fn cross(a: [f64; 3], b: [f64; 3]) -> [f64; 3] {
78    [
79        a[1] * b[2] - a[2] * b[1],
80        a[2] * b[0] - a[0] * b[2],
81        a[0] * b[1] - a[1] * b[0],
82    ]
83}
84
85fn dot(a: [f64; 3], b: [f64; 3]) -> f64 {
86    a[0] * b[0] + a[1] * b[1] + a[2] * b[2]
87}
88
89fn validate_particle_state(particle: &ChargedParticle, label: &str) -> FusionResult<()> {
90    if !particle.x_m.is_finite() || !particle.y_m.is_finite() || !particle.z_m.is_finite() {
91        return Err(FusionError::PhysicsViolation(format!(
92            "{label} position components must be finite"
93        )));
94    }
95    if !particle.vx_m_s.is_finite() || !particle.vy_m_s.is_finite() || !particle.vz_m_s.is_finite()
96    {
97        return Err(FusionError::PhysicsViolation(format!(
98            "{label} velocity components must be finite"
99        )));
100    }
101    if !particle.charge_c.is_finite() {
102        return Err(FusionError::PhysicsViolation(format!(
103            "{label}.charge_c must be finite"
104        )));
105    }
106    if !particle.mass_kg.is_finite() || particle.mass_kg <= 0.0 {
107        return Err(FusionError::PhysicsViolation(format!(
108            "{label}.mass_kg must be finite and > 0"
109        )));
110    }
111    if !particle.weight.is_finite() || particle.weight <= 0.0 {
112        return Err(FusionError::PhysicsViolation(format!(
113            "{label}.weight must be finite and > 0"
114        )));
115    }
116    Ok(())
117}
118
119fn validate_particle_projection_grid(grid: &Grid2D, label: &str) -> FusionResult<()> {
120    if grid.nr == 0 || grid.nz == 0 {
121        return Err(FusionError::PhysicsViolation(format!(
122            "{label} requires non-empty grid dimensions"
123        )));
124    }
125    if grid.r.len() != grid.nr || grid.z.len() != grid.nz {
126        return Err(FusionError::PhysicsViolation(format!(
127            "{label} grid axis length mismatch: r_len={}, nr={}, z_len={}, nz={}",
128            grid.r.len(),
129            grid.nr,
130            grid.z.len(),
131            grid.nz
132        )));
133    }
134    if grid.rr.dim() != (grid.nz, grid.nr) || grid.zz.dim() != (grid.nz, grid.nr) {
135        return Err(FusionError::PhysicsViolation(format!(
136            "{label} grid mesh shape mismatch: rr={:?}, zz={:?}, expected=({}, {})",
137            grid.rr.dim(),
138            grid.zz.dim(),
139            grid.nz,
140            grid.nr
141        )));
142    }
143    if grid.rr.iter().any(|v| !v.is_finite()) || grid.zz.iter().any(|v| !v.is_finite()) {
144        return Err(FusionError::PhysicsViolation(format!(
145            "{label} grid mesh coordinates must be finite"
146        )));
147    }
148    if !grid.dr.is_finite() || !grid.dz.is_finite() || grid.dr == 0.0 || grid.dz == 0.0 {
149        return Err(FusionError::PhysicsViolation(format!(
150            "{label} grid spacing must be finite and non-zero, got dr={}, dz={}",
151            grid.dr, grid.dz
152        )));
153    }
154    if grid.r.iter().any(|v| !v.is_finite()) || grid.z.iter().any(|v| !v.is_finite()) {
155        return Err(FusionError::PhysicsViolation(format!(
156            "{label} grid axes must be finite"
157        )));
158    }
159    Ok(())
160}
161
162fn nearest_index(axis: &Array1<f64>, value: f64, label: &str) -> FusionResult<usize> {
163    if axis.is_empty() {
164        return Err(FusionError::PhysicsViolation(format!(
165            "{label} axis must be non-empty"
166        )));
167    }
168    if !value.is_finite() {
169        return Err(FusionError::PhysicsViolation(format!(
170            "{label} lookup coordinate must be finite"
171        )));
172    }
173    let mut best_idx = 0usize;
174    let mut best_dist = f64::INFINITY;
175    for (idx, x) in axis.iter().copied().enumerate() {
176        if !x.is_finite() {
177            return Err(FusionError::PhysicsViolation(format!(
178                "{label} axis contains non-finite coordinate at index {idx}"
179            )));
180        }
181        let dist = (x - value).abs();
182        if dist < best_dist {
183            best_dist = dist;
184            best_idx = idx;
185        }
186    }
187    Ok(best_idx)
188}
189
190/// Create deterministic alpha test particles for hybrid kinetic-fluid overlays.
191pub fn seed_alpha_test_particles(
192    n_particles: usize,
193    major_radius_m: f64,
194    z_m: f64,
195    kinetic_energy_mev: f64,
196    pitch_cos: f64,
197    weight_per_particle: f64,
198) -> FusionResult<Vec<ChargedParticle>> {
199    if n_particles == 0 {
200        return Err(FusionError::PhysicsViolation(
201            "n_particles must be >= 1".to_string(),
202        ));
203    }
204    if !major_radius_m.is_finite() || major_radius_m <= 0.0 {
205        return Err(FusionError::PhysicsViolation(
206            "major_radius_m must be finite and > 0".to_string(),
207        ));
208    }
209    if !z_m.is_finite() {
210        return Err(FusionError::PhysicsViolation(
211            "z_m must be finite".to_string(),
212        ));
213    }
214    if !kinetic_energy_mev.is_finite() || kinetic_energy_mev <= 0.0 {
215        return Err(FusionError::PhysicsViolation(
216            "kinetic_energy_mev must be finite and > 0".to_string(),
217        ));
218    }
219    if !pitch_cos.is_finite() || !(-1.0..=1.0).contains(&pitch_cos) {
220        return Err(FusionError::PhysicsViolation(
221            "pitch_cos must be finite and in [-1, 1]".to_string(),
222        ));
223    }
224    if !weight_per_particle.is_finite() || weight_per_particle <= 0.0 {
225        return Err(FusionError::PhysicsViolation(
226            "weight_per_particle must be finite and > 0".to_string(),
227        ));
228    }
229
230    let r0 = major_radius_m;
231    let energy_j = kinetic_energy_mev * 1.0e6 * ELEMENTARY_CHARGE_C;
232    if !energy_j.is_finite() || energy_j <= 0.0 {
233        return Err(FusionError::PhysicsViolation(format!(
234            "seeded particle energy_j must be finite and > 0, got {energy_j}"
235        )));
236    }
237    let speed = (2.0 * energy_j / ALPHA_MASS_KG).sqrt();
238    if !speed.is_finite() || speed <= 0.0 {
239        return Err(FusionError::PhysicsViolation(format!(
240            "seeded particle speed must be finite and > 0, got {speed}"
241        )));
242    }
243    let pitch = pitch_cos;
244    let v_par = speed * pitch;
245    let perp_factor = (1.0 - pitch * pitch).max(0.0);
246    let v_perp = speed * perp_factor.sqrt();
247    if !v_par.is_finite() || !v_perp.is_finite() {
248        return Err(FusionError::PhysicsViolation(
249            "seeded particle velocity components became non-finite".to_string(),
250        ));
251    }
252    let weight = weight_per_particle;
253
254    let mut out = Vec::with_capacity(n_particles);
255    for i in 0..n_particles {
256        let phi = 2.0 * PI * (i as f64) / (n_particles as f64);
257        let x = r0 * phi.cos();
258        let y = r0 * phi.sin();
259        let ex = -phi.sin();
260        let ey = phi.cos();
261        out.push(ChargedParticle {
262            x_m: x,
263            y_m: y,
264            z_m,
265            vx_m_s: v_perp * ex,
266            vy_m_s: v_perp * ey,
267            vz_m_s: v_par,
268            charge_c: ALPHA_CHARGE_C,
269            mass_kg: ALPHA_MASS_KG,
270            weight,
271        });
272    }
273    Ok(out)
274}
275
276/// Summarize kinetic state of a particle population.
277pub fn summarize_particle_population(
278    particles: &[ChargedParticle],
279    runaway_threshold_mev: f64,
280) -> FusionResult<ParticlePopulationSummary> {
281    if !runaway_threshold_mev.is_finite() || runaway_threshold_mev < 0.0 {
282        return Err(FusionError::PhysicsViolation(
283            "runaway_threshold_mev must be finite and >= 0".to_string(),
284        ));
285    }
286    if particles.is_empty() {
287        return Ok(ParticlePopulationSummary {
288            count: 0,
289            mean_energy_mev: 0.0,
290            p95_energy_mev: 0.0,
291            max_energy_mev: 0.0,
292            runaway_fraction: 0.0,
293        });
294    }
295
296    let mut energies: Vec<f64> = Vec::with_capacity(particles.len());
297    for (idx, particle) in particles.iter().enumerate() {
298        validate_particle_state(particle, &format!("particle[{idx}]"))?;
299        let energy = particle.kinetic_energy_mev();
300        if !energy.is_finite() || energy < 0.0 {
301            return Err(FusionError::PhysicsViolation(format!(
302                "particle[{idx}] kinetic energy must be finite and >= 0, got {energy}"
303            )));
304        }
305        energies.push(energy);
306    }
307    energies.sort_by(f64::total_cmp);
308    let count = energies.len();
309    let mean_energy_mev = energies.iter().sum::<f64>() / (count as f64);
310    if !mean_energy_mev.is_finite() {
311        return Err(FusionError::PhysicsViolation(
312            "mean particle energy became non-finite".to_string(),
313        ));
314    }
315    let p95_idx = ((count - 1) as f64 * 0.95).round() as usize;
316    let p95_energy_mev = energies[p95_idx.min(count - 1)];
317    let max_energy_mev = energies[count - 1];
318    let n_runaway = energies
319        .iter()
320        .filter(|&&e| e >= runaway_threshold_mev)
321        .count();
322    let runaway_fraction = (n_runaway as f64) / (count as f64);
323    if !runaway_fraction.is_finite() {
324        return Err(FusionError::PhysicsViolation(
325            "runaway fraction became non-finite".to_string(),
326        ));
327    }
328
329    Ok(ParticlePopulationSummary {
330        count,
331        mean_energy_mev,
332        p95_energy_mev,
333        max_energy_mev,
334        runaway_fraction,
335    })
336}
337
338/// Estimate alpha-particle heating power density [W/m^3] on the R-Z grid.
339pub fn estimate_alpha_heating_profile(
340    particles: &[ChargedParticle],
341    grid: &Grid2D,
342    confinement_tau_s: f64,
343) -> FusionResult<Array2<f64>> {
344    if !confinement_tau_s.is_finite() || confinement_tau_s <= 0.0 {
345        return Err(FusionError::PhysicsViolation(
346            "confinement_tau_s must be finite and > 0".to_string(),
347        ));
348    }
349    validate_particle_projection_grid(grid, "alpha heating profile")?;
350    let mut heat = Array2::zeros((grid.nz, grid.nr));
351    if particles.is_empty() {
352        return Ok(heat);
353    }
354
355    let tau = confinement_tau_s;
356    let cell_volume = (grid.dr.abs() * grid.dz.abs() * 2.0 * PI).max(MIN_CELL_AREA_M2);
357    let r_min = grid.r[0].min(grid.r[grid.nr - 1]);
358    let r_max = grid.r[0].max(grid.r[grid.nr - 1]);
359    let z_min = grid.z[0].min(grid.z[grid.nz - 1]);
360    let z_max = grid.z[0].max(grid.z[grid.nz - 1]);
361
362    for (idx, particle) in particles.iter().enumerate() {
363        validate_particle_state(particle, &format!("particle[{idx}]"))?;
364        let r = particle.cylindrical_radius_m();
365        let z = particle.z_m;
366        if r < r_min || r > r_max || z < z_min || z > z_max {
367            continue;
368        }
369        let ir = nearest_index(&grid.r, r, "alpha heating R-axis")?;
370        let iz = nearest_index(&grid.z, z, "alpha heating Z-axis")?;
371        let p_w = particle.kinetic_energy_j() * particle.weight / tau;
372        if !p_w.is_finite() {
373            return Err(FusionError::PhysicsViolation(format!(
374                "particle[{idx}] deposited heating power became non-finite"
375            )));
376        }
377        let local_volume = (cell_volume * r.max(MIN_RADIUS_M)).max(MIN_CELL_AREA_M2);
378        let contribution = p_w / local_volume;
379        if !contribution.is_finite() {
380            return Err(FusionError::PhysicsViolation(format!(
381                "particle[{idx}] heating contribution became non-finite"
382            )));
383        }
384        heat[[iz, ir]] += contribution;
385    }
386    if heat.iter().any(|v| !v.is_finite()) {
387        return Err(FusionError::PhysicsViolation(
388            "alpha heating profile contains non-finite values".to_string(),
389        ));
390    }
391    Ok(heat)
392}
393
394/// Advance one particle state using the Boris push.
395pub fn boris_push_step(
396    particle: &mut ChargedParticle,
397    electric_v_m: [f64; 3],
398    magnetic_t: [f64; 3],
399    dt_s: f64,
400) -> FusionResult<()> {
401    validate_particle_state(particle, "particle")?;
402    if electric_v_m.iter().any(|v| !v.is_finite()) || magnetic_t.iter().any(|v| !v.is_finite()) {
403        return Err(FusionError::PhysicsViolation(
404            "electric_v_m and magnetic_t must be finite vectors".to_string(),
405        ));
406    }
407    if !dt_s.is_finite() || dt_s <= 0.0 {
408        return Err(FusionError::PhysicsViolation(
409            "dt_s must be finite and > 0".to_string(),
410        ));
411    }
412
413    let qmdt2 = particle.charge_c * dt_s / (2.0 * particle.mass_kg);
414    if !qmdt2.is_finite() {
415        return Err(FusionError::PhysicsViolation(
416            "boris qmdt2 became non-finite".to_string(),
417        ));
418    }
419    let v_minus = [
420        particle.vx_m_s + qmdt2 * electric_v_m[0],
421        particle.vy_m_s + qmdt2 * electric_v_m[1],
422        particle.vz_m_s + qmdt2 * electric_v_m[2],
423    ];
424
425    let t = [
426        qmdt2 * magnetic_t[0],
427        qmdt2 * magnetic_t[1],
428        qmdt2 * magnetic_t[2],
429    ];
430    let t2 = dot(t, t);
431    let s = [
432        (2.0 * t[0]) / (1.0 + t2),
433        (2.0 * t[1]) / (1.0 + t2),
434        (2.0 * t[2]) / (1.0 + t2),
435    ];
436
437    let v_prime = {
438        let c = cross(v_minus, t);
439        [v_minus[0] + c[0], v_minus[1] + c[1], v_minus[2] + c[2]]
440    };
441    let v_plus = {
442        let c = cross(v_prime, s);
443        [v_minus[0] + c[0], v_minus[1] + c[1], v_minus[2] + c[2]]
444    };
445
446    let vx_new = v_plus[0] + qmdt2 * electric_v_m[0];
447    let vy_new = v_plus[1] + qmdt2 * electric_v_m[1];
448    let vz_new = v_plus[2] + qmdt2 * electric_v_m[2];
449    if !vx_new.is_finite() || !vy_new.is_finite() || !vz_new.is_finite() {
450        return Err(FusionError::PhysicsViolation(
451            "boris velocity update became non-finite".to_string(),
452        ));
453    }
454
455    particle.vx_m_s = vx_new;
456    particle.vy_m_s = vy_new;
457    particle.vz_m_s = vz_new;
458    particle.x_m += vx_new * dt_s;
459    particle.y_m += vy_new * dt_s;
460    particle.z_m += vz_new * dt_s;
461    if !particle.x_m.is_finite() || !particle.y_m.is_finite() || !particle.z_m.is_finite() {
462        return Err(FusionError::PhysicsViolation(
463            "boris position update became non-finite".to_string(),
464        ));
465    }
466    Ok(())
467}
468
469/// Advance a particle set for a fixed number of Boris steps.
470pub fn advance_particles_boris(
471    particles: &mut [ChargedParticle],
472    electric_v_m: [f64; 3],
473    magnetic_t: [f64; 3],
474    dt_s: f64,
475    steps: usize,
476) -> FusionResult<()> {
477    if particles.is_empty() {
478        return Err(FusionError::PhysicsViolation(
479            "particles must be non-empty".to_string(),
480        ));
481    }
482    if steps == 0 {
483        return Err(FusionError::PhysicsViolation(
484            "steps must be >= 1".to_string(),
485        ));
486    }
487    if !dt_s.is_finite() || dt_s <= 0.0 {
488        return Err(FusionError::PhysicsViolation(
489            "dt_s must be finite and > 0".to_string(),
490        ));
491    }
492    for _ in 0..steps {
493        for (idx, particle) in particles.iter_mut().enumerate() {
494            boris_push_step(particle, electric_v_m, magnetic_t, dt_s).map_err(|err| match err {
495                FusionError::PhysicsViolation(msg) => FusionError::PhysicsViolation(format!(
496                    "particle[{idx}] boris push failed: {msg}"
497                )),
498                other => other,
499            })?;
500        }
501    }
502    Ok(())
503}
504
505/// Deposit particle toroidal current density J_phi on the GS R-Z grid.
506pub fn deposit_toroidal_current_density(
507    particles: &[ChargedParticle],
508    grid: &Grid2D,
509) -> FusionResult<Array2<f64>> {
510    validate_particle_projection_grid(grid, "particle current deposition")?;
511    if particles.is_empty() {
512        return Err(FusionError::PhysicsViolation(
513            "particles must be non-empty".to_string(),
514        ));
515    }
516    let mut j_phi: Array2<f64> = Array2::zeros((grid.nz, grid.nr));
517    let area = (grid.dr.abs() * grid.dz.abs()).max(MIN_CELL_AREA_M2);
518    let r_min = grid.r[0].min(grid.r[grid.nr - 1]);
519    let r_max = grid.r[0].max(grid.r[grid.nr - 1]);
520    let z_min = grid.z[0].min(grid.z[grid.nz - 1]);
521    let z_max = grid.z[0].max(grid.z[grid.nz - 1]);
522
523    for (idx, particle) in particles.iter().enumerate() {
524        validate_particle_state(particle, &format!("particle[{idx}]"))?;
525        let r = particle.cylindrical_radius_m();
526        let z = particle.z_m;
527        if r < r_min || r > r_max || z < z_min || z > z_max {
528            continue;
529        }
530
531        let ir = nearest_index(&grid.r, r, "particle current R-axis")?;
532        let iz = nearest_index(&grid.z, z, "particle current Z-axis")?;
533        let v_phi = particle.toroidal_velocity_m_s();
534        let j_contrib = particle.charge_c * particle.weight * v_phi / area;
535        if !j_contrib.is_finite() {
536            return Err(FusionError::PhysicsViolation(format!(
537                "particle[{idx}] toroidal current contribution became non-finite"
538            )));
539        }
540        j_phi[[iz, ir]] += j_contrib;
541    }
542
543    if j_phi.iter().any(|v| !v.is_finite()) {
544        return Err(FusionError::PhysicsViolation(
545            "deposited particle current contains non-finite values".to_string(),
546        ));
547    }
548    Ok(j_phi)
549}
550
551/// Blend fluid and particle current maps and renormalize integral to `i_target`.
552pub fn blend_particle_current(
553    fluid_j_phi: &Array2<f64>,
554    particle_j_phi: &Array2<f64>,
555    grid: &Grid2D,
556    i_target: f64,
557    particle_coupling: f64,
558) -> FusionResult<Array2<f64>> {
559    validate_particle_projection_grid(grid, "particle current blend")?;
560    let expected_shape = (grid.nz, grid.nr);
561    if fluid_j_phi.dim() != expected_shape || particle_j_phi.dim() != expected_shape {
562        return Err(FusionError::PhysicsViolation(format!(
563            "Particle-current blend shape mismatch: expected {:?}, fluid {:?}, particle {:?}",
564            expected_shape,
565            fluid_j_phi.dim(),
566            particle_j_phi.dim(),
567        )));
568    }
569    if !particle_coupling.is_finite() || !(0.0..=1.0).contains(&particle_coupling) {
570        return Err(FusionError::PhysicsViolation(
571            "particle_coupling must be finite and in [0, 1]".to_string(),
572        ));
573    }
574    if !i_target.is_finite() {
575        return Err(FusionError::PhysicsViolation(
576            "i_target must be finite".to_string(),
577        ));
578    }
579    if fluid_j_phi.iter().any(|v| !v.is_finite()) || particle_j_phi.iter().any(|v| !v.is_finite()) {
580        return Err(FusionError::PhysicsViolation(
581            "fluid_j_phi and particle_j_phi must be finite".to_string(),
582        ));
583    }
584
585    let coupling = particle_coupling;
586    let fluid_weight = 1.0 - coupling;
587    let mut combined = Array2::zeros(expected_shape);
588
589    for iz in 0..grid.nz {
590        for ir in 0..grid.nr {
591            combined[[iz, ir]] =
592                fluid_weight * fluid_j_phi[[iz, ir]] + coupling * particle_j_phi[[iz, ir]];
593        }
594    }
595
596    let i_current = combined.iter().sum::<f64>() * grid.dr * grid.dz;
597    if !i_current.is_finite() {
598        return Err(FusionError::PhysicsViolation(
599            "blended current integral became non-finite".to_string(),
600        ));
601    }
602    if i_current.abs() > MIN_CURRENT_INTEGRAL {
603        let scale = i_target / i_current;
604        if !scale.is_finite() {
605            return Err(FusionError::PhysicsViolation(
606                "blended current scale became non-finite".to_string(),
607            ));
608        }
609        combined.mapv_inplace(|v| v * scale);
610    } else if i_target.abs() > MIN_CURRENT_INTEGRAL {
611        return Err(FusionError::PhysicsViolation(format!(
612            "cannot renormalize blended current: |i_current| <= {MIN_CURRENT_INTEGRAL} while |i_target| > {MIN_CURRENT_INTEGRAL}"
613        )));
614    } else {
615        combined.fill(0.0);
616    }
617    if combined.iter().any(|v| !v.is_finite()) {
618        return Err(FusionError::PhysicsViolation(
619            "blended current map contains non-finite values".to_string(),
620        ));
621    }
622
623    Ok(combined)
624}
625
626// ═══════════════════════════════════════════════════════════════════════
627// Coulomb Collision Operator (Fokker-Planck Monte Carlo)
628// ═══════════════════════════════════════════════════════════════════════
629
630/// Parameters describing the background plasma for Coulomb collisions.
631#[derive(Debug, Clone, Copy, PartialEq)]
632pub struct CoulombCollisionParams {
633    /// Electron density [m^-3].
634    pub n_e: f64,
635    /// Electron temperature [keV].
636    pub t_e_kev: f64,
637    /// Ion temperature [keV].
638    pub t_i_kev: f64,
639    /// Ion mass number (e.g. 2 for deuterium).
640    pub a_i: f64,
641    /// Ion charge number.
642    pub z_i: f64,
643    /// Effective charge Z_eff.
644    pub z_eff: f64,
645}
646
647fn validate_collision_params(p: &CoulombCollisionParams) -> FusionResult<()> {
648    if !p.n_e.is_finite() || p.n_e <= 0.0 {
649        return Err(FusionError::PhysicsViolation(
650            "n_e must be finite and > 0".into(),
651        ));
652    }
653    if !p.t_e_kev.is_finite() || p.t_e_kev <= 0.0 {
654        return Err(FusionError::PhysicsViolation(
655            "t_e_kev must be finite and > 0".into(),
656        ));
657    }
658    if !p.t_i_kev.is_finite() || p.t_i_kev <= 0.0 {
659        return Err(FusionError::PhysicsViolation(
660            "t_i_kev must be finite and > 0".into(),
661        ));
662    }
663    if !p.a_i.is_finite() || p.a_i <= 0.0 {
664        return Err(FusionError::PhysicsViolation(
665            "a_i must be finite and > 0".into(),
666        ));
667    }
668    if !p.z_i.is_finite() || p.z_i <= 0.0 {
669        return Err(FusionError::PhysicsViolation(
670            "z_i must be finite and > 0".into(),
671        ));
672    }
673    if !p.z_eff.is_finite() || p.z_eff < 1.0 {
674        return Err(FusionError::PhysicsViolation(
675            "z_eff must be finite and >= 1".into(),
676        ));
677    }
678    Ok(())
679}
680
681/// Coulomb logarithm via NRL formula, clamped to [5, 30].
682pub fn coulomb_logarithm(n_e_m3: f64, t_e_kev: f64) -> FusionResult<f64> {
683    if !n_e_m3.is_finite() || n_e_m3 <= 0.0 {
684        return Err(FusionError::PhysicsViolation(
685            "n_e must be finite and > 0".into(),
686        ));
687    }
688    if !t_e_kev.is_finite() || t_e_kev <= 0.0 {
689        return Err(FusionError::PhysicsViolation(
690            "t_e_kev must be finite and > 0".into(),
691        ));
692    }
693    // NRL Formulary: ln Λ ≈ 24 − ln(sqrt(n_e) / T_e)  [T_e in eV]
694    let t_e_ev = t_e_kev * 1000.0;
695    let ln_lambda = 24.0 - (n_e_m3.sqrt() / t_e_ev).ln();
696    Ok(ln_lambda.clamp(5.0, 30.0))
697}
698
699/// Spitzer slowing-down time [s] for a test particle on field electrons.
700///
701/// τ_s = 3(2π)^{3/2} ε₀² m_a T_e^{3/2} / (n_e Z_a² e⁴ m_e^{1/2} ln Λ)
702pub fn spitzer_slowing_down_time(
703    mass_kg: f64,
704    charge_number: f64,
705    n_e_m3: f64,
706    t_e_kev: f64,
707    ln_lambda: f64,
708) -> FusionResult<f64> {
709    if !mass_kg.is_finite() || mass_kg <= 0.0 {
710        return Err(FusionError::PhysicsViolation(
711            "mass_kg must be finite and > 0".into(),
712        ));
713    }
714    if !charge_number.is_finite() || charge_number <= 0.0 {
715        return Err(FusionError::PhysicsViolation(
716            "charge_number must be > 0".into(),
717        ));
718    }
719    if !n_e_m3.is_finite() || n_e_m3 <= 0.0 {
720        return Err(FusionError::PhysicsViolation("n_e must be > 0".into()));
721    }
722    if !t_e_kev.is_finite() || t_e_kev <= 0.0 {
723        return Err(FusionError::PhysicsViolation("t_e_kev must be > 0".into()));
724    }
725    if !ln_lambda.is_finite() || ln_lambda <= 0.0 {
726        return Err(FusionError::PhysicsViolation(
727            "ln_lambda must be > 0".into(),
728        ));
729    }
730    let t_e_j = t_e_kev * BOLTZMANN_J_PER_KEV;
731    let e = ELEMENTARY_CHARGE_C;
732    let numerator = 3.0
733        * (2.0 * PI).powf(1.5)
734        * VACUUM_PERMITTIVITY
735        * VACUUM_PERMITTIVITY
736        * mass_kg
737        * t_e_j.powf(1.5);
738    let denominator =
739        n_e_m3 * (charge_number * e).powi(2) * e * e * ELECTRON_MASS_KG.sqrt() * ln_lambda;
740    let tau = numerator / denominator;
741    if !tau.is_finite() || tau <= 0.0 {
742        return Err(FusionError::PhysicsViolation(format!(
743            "spitzer time became non-physical: {tau}"
744        )));
745    }
746    Ok(tau)
747}
748
749/// Critical velocity where electron drag equals ion drag [m/s].
750pub fn critical_velocity(t_e_kev: f64, a_i: f64, z_i: f64, z_eff: f64) -> FusionResult<f64> {
751    if !t_e_kev.is_finite() || t_e_kev <= 0.0 {
752        return Err(FusionError::PhysicsViolation("t_e_kev must be > 0".into()));
753    }
754    if !a_i.is_finite() || a_i <= 0.0 {
755        return Err(FusionError::PhysicsViolation("a_i must be > 0".into()));
756    }
757    if !z_i.is_finite() || z_i <= 0.0 {
758        return Err(FusionError::PhysicsViolation("z_i must be > 0".into()));
759    }
760    if !z_eff.is_finite() || z_eff < 1.0 {
761        return Err(FusionError::PhysicsViolation("z_eff must be >= 1".into()));
762    }
763    let t_e_j = t_e_kev * BOLTZMANN_J_PER_KEV;
764    let v_te = (2.0 * t_e_j / ELECTRON_MASS_KG).sqrt();
765    // v_c = v_te * (3 sqrt(π) m_e / (4 m_i))^{1/3} * Z_eff^{1/3}
766    let mass_ratio = ELECTRON_MASS_KG / (a_i * PROTON_MASS_KG);
767    let factor = (0.75 * PI.sqrt() * mass_ratio).powf(1.0 / 3.0);
768    let v_c = v_te * factor * z_eff.powf(1.0 / 3.0);
769    if !v_c.is_finite() || v_c <= 0.0 {
770        return Err(FusionError::PhysicsViolation(format!(
771            "critical velocity non-physical: {v_c}"
772        )));
773    }
774    Ok(v_c)
775}
776
777/// Collision frequencies (ν_slow, ν_defl, ν_energy) for a test particle at given speed.
778pub fn collision_frequencies(
779    speed: f64,
780    params: &CoulombCollisionParams,
781    ln_lambda: f64,
782    tau_s: f64,
783    v_c: f64,
784) -> FusionResult<(f64, f64, f64)> {
785    if !speed.is_finite() || speed < 0.0 {
786        return Err(FusionError::PhysicsViolation(
787            "speed must be finite and >= 0".into(),
788        ));
789    }
790    if !tau_s.is_finite() || tau_s <= 0.0 {
791        return Err(FusionError::PhysicsViolation("tau_s must be > 0".into()));
792    }
793    if !v_c.is_finite() || v_c <= 0.0 {
794        return Err(FusionError::PhysicsViolation("v_c must be > 0".into()));
795    }
796    let _ = ln_lambda; // used indirectly via tau_s
797    let _ = params;
798
799    let v_safe = speed.max(1e-6);
800    let x3 = (v_safe / v_c).powi(3);
801
802    // Slowing-down frequency: ν_s = (1 + x³) / τ_s  where x = v/v_c
803    let nu_slow = (1.0 + x3) / tau_s;
804
805    // Deflection (pitch-angle scattering): ν_d ≈ (1 + Z_eff/2) / (τ_s · x³)
806    let nu_defl = if x3 > 1e-30 {
807        1.0 / (tau_s * x3)
808    } else {
809        1.0 / (tau_s * 1e-30)
810    };
811
812    // Energy diffusion: ν_ε ≈ 2 T_e / (m v² τ_s)
813    let nu_energy = nu_slow * 0.5; // simplified
814
815    Ok((nu_slow, nu_defl, nu_energy))
816}
817
818/// XORshift64 PRNG: returns next state and a uniform f64 in [0, 1).
819fn xorshift_uniform(state: &mut u64) -> f64 {
820    let mut s = *state;
821    s ^= s << 13;
822    s ^= s >> 7;
823    s ^= s << 17;
824    *state = s;
825    (s as f64) / (u64::MAX as f64)
826}
827
828/// Box-Muller transform producing one standard normal variate from xorshift.
829fn xorshift_normal(state: &mut u64) -> f64 {
830    let u1 = xorshift_uniform(state).max(1e-300);
831    let u2 = xorshift_uniform(state);
832    (-2.0 * u1.ln()).sqrt() * (2.0 * PI * u2).cos()
833}
834
835/// Find an orthonormal basis (e1, e2) perpendicular to unit vector `v_hat`.
836fn perpendicular_basis(v_hat: [f64; 3]) -> ([f64; 3], [f64; 3]) {
837    // Choose the axis most perpendicular to v_hat
838    let abs_x = v_hat[0].abs();
839    let abs_y = v_hat[1].abs();
840    let abs_z = v_hat[2].abs();
841    let seed = if abs_x <= abs_y && abs_x <= abs_z {
842        [1.0, 0.0, 0.0]
843    } else if abs_y <= abs_z {
844        [0.0, 1.0, 0.0]
845    } else {
846        [0.0, 0.0, 1.0]
847    };
848
849    // e1 = normalize(seed × v_hat)
850    let raw = cross(seed, v_hat);
851    let norm = dot(raw, raw).sqrt().max(1e-30);
852    let e1 = [raw[0] / norm, raw[1] / norm, raw[2] / norm];
853
854    // e2 = v_hat × e1
855    let e2 = cross(v_hat, e1);
856    (e1, e2)
857}
858
859/// Apply Coulomb collision Monte Carlo kick to a single particle.
860///
861/// Uses Langevin approach: drag (slowing-down) + stochastic pitch-angle
862/// scattering and energy diffusion over timestep dt_s.
863pub fn collision_step(
864    particle: &mut ChargedParticle,
865    params: &CoulombCollisionParams,
866    dt_s: f64,
867    rng_state: &mut u64,
868) -> FusionResult<()> {
869    validate_particle_state(particle, "collision particle")?;
870    validate_collision_params(params)?;
871    if !dt_s.is_finite() || dt_s <= 0.0 {
872        return Err(FusionError::PhysicsViolation("dt_s must be > 0".into()));
873    }
874
875    let v = [particle.vx_m_s, particle.vy_m_s, particle.vz_m_s];
876    let speed = dot(v, v).sqrt();
877    if speed < 1e-10 {
878        return Ok(()); // particle at rest — no collision
879    }
880
881    let ln_lam = coulomb_logarithm(params.n_e, params.t_e_kev)?;
882    let za = (particle.charge_c / ELEMENTARY_CHARGE_C).abs();
883    let tau_s =
884        spitzer_slowing_down_time(particle.mass_kg, za, params.n_e, params.t_e_kev, ln_lam)?;
885    let v_c = critical_velocity(params.t_e_kev, params.a_i, params.z_i, params.z_eff)?;
886    let (nu_slow, nu_defl, _nu_energy) = collision_frequencies(speed, params, ln_lam, tau_s, v_c)?;
887
888    // Unit velocity direction
889    let v_hat = [v[0] / speed, v[1] / speed, v[2] / speed];
890
891    // 1. Slowing-down drag: Δv_∥ = -ν_s · v · dt
892    let dv_par = -nu_slow * speed * dt_s;
893
894    // 2. Pitch-angle scattering: random perpendicular kick
895    //    σ_perp = v · sqrt(ν_d · dt)
896    let sigma_perp = speed * (nu_defl * dt_s).sqrt();
897    let kick1 = xorshift_normal(rng_state) * sigma_perp;
898    let kick2 = xorshift_normal(rng_state) * sigma_perp;
899
900    let (e1, e2) = perpendicular_basis(v_hat);
901
902    // New velocity = (v + dv_par) v_hat + kick1 e1 + kick2 e2
903    let new_speed_par = speed + dv_par;
904    // Ensure speed doesn't go negative (particle thermalized)
905    let new_speed_par = new_speed_par.max(0.0);
906
907    particle.vx_m_s = new_speed_par * v_hat[0] + kick1 * e1[0] + kick2 * e2[0];
908    particle.vy_m_s = new_speed_par * v_hat[1] + kick1 * e1[1] + kick2 * e2[1];
909    particle.vz_m_s = new_speed_par * v_hat[2] + kick1 * e1[2] + kick2 * e2[2];
910
911    if !particle.vx_m_s.is_finite() || !particle.vy_m_s.is_finite() || !particle.vz_m_s.is_finite()
912    {
913        return Err(FusionError::PhysicsViolation(
914            "collision step produced non-finite velocity".into(),
915        ));
916    }
917    Ok(())
918}
919
920/// Apply Coulomb collisions to all particles in a batch.
921pub fn apply_coulomb_collisions(
922    particles: &mut [ChargedParticle],
923    params: &CoulombCollisionParams,
924    dt_s: f64,
925    seed: u64,
926) -> FusionResult<()> {
927    validate_collision_params(params)?;
928    if particles.is_empty() {
929        return Ok(());
930    }
931    if !dt_s.is_finite() || dt_s <= 0.0 {
932        return Err(FusionError::PhysicsViolation("dt_s must be > 0".into()));
933    }
934    if seed == 0 {
935        return Err(FusionError::PhysicsViolation(
936            "seed must be != 0 for xorshift".into(),
937        ));
938    }
939
940    for (idx, particle) in particles.iter_mut().enumerate() {
941        // Each particle gets a deterministic RNG state derived from seed + index
942        let mut rng = seed
943            .wrapping_add(idx as u64)
944            .wrapping_mul(6364136223846793005)
945            .wrapping_add(1);
946        if rng == 0 {
947            rng = 1;
948        }
949        collision_step(particle, params, dt_s, &mut rng).map_err(|e| match e {
950            FusionError::PhysicsViolation(msg) => {
951                FusionError::PhysicsViolation(format!("particle[{idx}] collision failed: {msg}"))
952            }
953            other => other,
954        })?;
955    }
956    Ok(())
957}
958
959#[cfg(test)]
960mod tests {
961    use super::*;
962
963    #[test]
964    fn test_boris_push_preserves_speed_without_electric_field() {
965        let particle = ChargedParticle {
966            x_m: 2.0,
967            y_m: 0.0,
968            z_m: 0.0,
969            vx_m_s: 80_000.0,
970            vy_m_s: 10_000.0,
971            vz_m_s: 5_000.0,
972            charge_c: 1.602_176_634e-19,
973            mass_kg: 1.672_621_923_69e-27,
974            weight: 1.0,
975        };
976        let speed_0 =
977            (particle.vx_m_s.powi(2) + particle.vy_m_s.powi(2) + particle.vz_m_s.powi(2)).sqrt();
978        let mut particles = vec![particle];
979        advance_particles_boris(&mut particles, [0.0, 0.0, 0.0], [0.0, 0.0, 2.5], 5e-10, 600)
980            .expect("valid boris advance should succeed");
981        let updated = particles[0];
982        let speed_1 =
983            (updated.vx_m_s.powi(2) + updated.vy_m_s.powi(2) + updated.vz_m_s.powi(2)).sqrt();
984        let rel = (speed_1 - speed_0).abs() / speed_0;
985        assert!(rel < 1e-10, "Speed drift too high for Boris push: {rel}");
986    }
987
988    #[test]
989    fn test_toroidal_current_deposition_is_nonzero() {
990        let grid = Grid2D::new(17, 17, 1.0, 9.0, -4.0, 4.0);
991        let particles = vec![
992            ChargedParticle {
993                x_m: 5.0,
994                y_m: 0.0,
995                z_m: 0.0,
996                vx_m_s: 0.0,
997                vy_m_s: 100_000.0,
998                vz_m_s: 0.0,
999                charge_c: 1.602_176_634e-19,
1000                mass_kg: 1.672_621_923_69e-27,
1001                weight: 2.0e16,
1002            },
1003            ChargedParticle {
1004                x_m: 5.0,
1005                y_m: 0.2,
1006                z_m: 0.3,
1007                vx_m_s: 10_000.0,
1008                vy_m_s: 90_000.0,
1009                vz_m_s: -2_000.0,
1010                charge_c: 1.602_176_634e-19,
1011                mass_kg: 1.672_621_923_69e-27,
1012                weight: 1.5e16,
1013            },
1014        ];
1015
1016        let j = deposit_toroidal_current_density(&particles, &grid)
1017            .expect("valid particles/grid should deposit current");
1018        let sum_abs = j.iter().map(|v| v.abs()).sum::<f64>();
1019        assert!(
1020            sum_abs > 0.0,
1021            "Expected non-zero toroidal current deposition"
1022        );
1023    }
1024
1025    #[test]
1026    fn test_blend_particle_current_renormalizes_target_current() {
1027        let grid = Grid2D::new(8, 8, 1.0, 5.0, -2.0, 2.0);
1028        let fluid = Array2::from_elem((8, 8), 2.0);
1029        let particle = Array2::from_elem((8, 8), 6.0);
1030        let i_target = 15.0e6;
1031        let blended = blend_particle_current(&fluid, &particle, &grid, i_target, 0.25).unwrap();
1032        let i_actual = blended.iter().sum::<f64>() * grid.dr * grid.dz;
1033        let rel = ((i_actual - i_target) / i_target).abs();
1034        assert!(
1035            rel < 1e-12,
1036            "Blended current should match target after renormalization"
1037        );
1038    }
1039
1040    #[test]
1041    fn test_blend_particle_current_shape_mismatch_errors() {
1042        let grid = Grid2D::new(8, 8, 1.0, 5.0, -2.0, 2.0);
1043        let fluid = Array2::zeros((8, 8));
1044        let particle_bad = Array2::zeros((7, 8));
1045        let err = blend_particle_current(&fluid, &particle_bad, &grid, 1.0, 0.5).unwrap_err();
1046        match err {
1047            FusionError::PhysicsViolation(msg) => {
1048                assert!(msg.contains("shape mismatch"));
1049            }
1050            other => panic!("Unexpected error: {other:?}"),
1051        }
1052    }
1053
1054    #[test]
1055    fn test_blend_particle_current_rejects_invalid_coupling_and_target() {
1056        let grid = Grid2D::new(8, 8, 1.0, 5.0, -2.0, 2.0);
1057        let fluid = Array2::from_elem((8, 8), 2.0);
1058        let particle = Array2::from_elem((8, 8), 6.0);
1059        for bad_coupling in [f64::NAN, -0.1, 1.1] {
1060            let err = blend_particle_current(&fluid, &particle, &grid, 1.0, bad_coupling)
1061                .expect_err("invalid coupling must error");
1062            match err {
1063                FusionError::PhysicsViolation(msg) => {
1064                    assert!(msg.contains("particle_coupling"));
1065                }
1066                other => panic!("Unexpected error: {other:?}"),
1067            }
1068        }
1069        let err = blend_particle_current(&fluid, &particle, &grid, f64::INFINITY, 0.5)
1070            .expect_err("non-finite target current must error");
1071        match err {
1072            FusionError::PhysicsViolation(msg) => {
1073                assert!(msg.contains("i_target"));
1074            }
1075            other => panic!("Unexpected error: {other:?}"),
1076        }
1077    }
1078
1079    #[test]
1080    fn test_seed_alpha_particles_matches_requested_energy_band() {
1081        let particles =
1082            seed_alpha_test_particles(24, 6.2, 0.0, 3.5, 0.6, 5.0e12).expect("valid seeds");
1083        assert_eq!(particles.len(), 24);
1084        let summary = summarize_particle_population(&particles, 2.0).expect("valid threshold");
1085        assert!(summary.mean_energy_mev > 3.0);
1086        assert!(summary.max_energy_mev < 4.2);
1087        assert!(summary.runaway_fraction > 0.9);
1088    }
1089
1090    #[test]
1091    fn test_alpha_heating_profile_is_positive_when_particles_in_domain() {
1092        let grid = Grid2D::new(33, 33, 3.0, 9.0, -2.5, 2.5);
1093        let particles =
1094            seed_alpha_test_particles(16, 6.0, 0.1, 3.5, 0.4, 1.0e13).expect("valid seeds");
1095        let heat =
1096            estimate_alpha_heating_profile(&particles, &grid, 0.25).expect("valid confinement");
1097        let total = heat.iter().sum::<f64>();
1098        assert!(total > 0.0, "Expected positive deposited alpha heating");
1099        assert!(!heat.iter().any(|v| !v.is_finite()));
1100    }
1101
1102    #[test]
1103    fn test_alpha_heating_profile_supports_descending_axes() {
1104        let grid = Grid2D::new(33, 33, 9.0, 3.0, 2.5, -2.5);
1105        let particles =
1106            seed_alpha_test_particles(16, 6.0, 0.1, 3.5, 0.4, 1.0e13).expect("valid seeds");
1107        let heat =
1108            estimate_alpha_heating_profile(&particles, &grid, 0.25).expect("valid confinement");
1109        let total = heat.iter().sum::<f64>();
1110        assert!(
1111            total > 0.0,
1112            "Expected positive deposited alpha heating on descending axes"
1113        );
1114    }
1115
1116    #[test]
1117    fn test_alpha_heating_profile_rejects_invalid_confinement_time() {
1118        let grid = Grid2D::new(17, 17, 3.0, 9.0, -1.5, 1.5);
1119        let particles =
1120            seed_alpha_test_particles(8, 6.0, 0.0, 3.5, 0.2, 1.0e12).expect("valid seeds");
1121        for bad_tau in [0.0, f64::NAN] {
1122            let err = estimate_alpha_heating_profile(&particles, &grid, bad_tau)
1123                .expect_err("invalid confinement time must error");
1124            match err {
1125                FusionError::PhysicsViolation(msg) => {
1126                    assert!(msg.contains("confinement_tau_s"));
1127                }
1128                other => panic!("Unexpected error: {other:?}"),
1129            }
1130        }
1131    }
1132
1133    #[test]
1134    fn test_blend_particle_current_rejects_non_finite_maps() {
1135        let grid = Grid2D::new(8, 8, 1.0, 5.0, -2.0, 2.0);
1136        let mut fluid = Array2::from_elem((8, 8), 2.0);
1137        let particle = Array2::from_elem((8, 8), 6.0);
1138        fluid[[0, 0]] = f64::NAN;
1139        let err = blend_particle_current(&fluid, &particle, &grid, 1.0, 0.5)
1140            .expect_err("non-finite current map must fail");
1141        match err {
1142            FusionError::PhysicsViolation(msg) => {
1143                assert!(msg.contains("must be finite"));
1144            }
1145            other => panic!("Unexpected error: {other:?}"),
1146        }
1147    }
1148
1149    #[test]
1150    fn test_blend_particle_current_rejects_zero_integral_nonzero_target() {
1151        let grid = Grid2D::new(8, 8, 1.0, 5.0, -2.0, 2.0);
1152        let fluid = Array2::zeros((8, 8));
1153        let particle = Array2::zeros((8, 8));
1154        let err = blend_particle_current(&fluid, &particle, &grid, 1.0e6, 0.5)
1155            .expect_err("non-zero target with zero blended current must fail");
1156        match err {
1157            FusionError::PhysicsViolation(msg) => {
1158                assert!(msg.contains("cannot renormalize blended current"));
1159            }
1160            other => panic!("Unexpected error: {other:?}"),
1161        }
1162    }
1163
1164    #[test]
1165    fn test_blend_particle_current_rejects_invalid_grid_spacing() {
1166        let mut grid = Grid2D::new(8, 8, 1.0, 5.0, -2.0, 2.0);
1167        grid.dr = 0.0;
1168        let fluid = Array2::from_elem((8, 8), 2.0);
1169        let particle = Array2::from_elem((8, 8), 6.0);
1170        let err = blend_particle_current(&fluid, &particle, &grid, 1.0, 0.5)
1171            .expect_err("invalid grid spacing must fail");
1172        match err {
1173            FusionError::PhysicsViolation(msg) => {
1174                assert!(msg.contains("grid spacing"));
1175            }
1176            other => panic!("Unexpected error: {other:?}"),
1177        }
1178    }
1179
1180    #[test]
1181    fn test_seed_alpha_particles_rejects_invalid_parameters() {
1182        let bad = [
1183            seed_alpha_test_particles(0, 6.2, 0.0, 3.5, 0.6, 1.0),
1184            seed_alpha_test_particles(8, 0.0, 0.0, 3.5, 0.6, 1.0),
1185            seed_alpha_test_particles(8, 6.2, f64::NAN, 3.5, 0.6, 1.0),
1186            seed_alpha_test_particles(8, 6.2, 0.0, -1.0, 0.6, 1.0),
1187            seed_alpha_test_particles(8, 6.2, 0.0, 3.5, 1.2, 1.0),
1188            seed_alpha_test_particles(8, 6.2, 0.0, 3.5, 0.6, 0.0),
1189        ];
1190        for candidate in bad {
1191            assert!(
1192                candidate.is_err(),
1193                "Expected invalid seed parameters to return an error"
1194            );
1195        }
1196    }
1197
1198    #[test]
1199    fn test_seed_alpha_particles_rejects_non_finite_kinematics() {
1200        let err = seed_alpha_test_particles(4, 6.2, 0.0, f64::MAX, 0.6, 1.0)
1201            .expect_err("overflowing kinetic seed should fail");
1202        match err {
1203            FusionError::PhysicsViolation(msg) => {
1204                assert!(msg.contains("energy_j") || msg.contains("speed"));
1205            }
1206            other => panic!("Unexpected error: {other:?}"),
1207        }
1208    }
1209
1210    #[test]
1211    fn test_population_summary_empty_is_zero() {
1212        let summary = summarize_particle_population(&[], 1.0).expect("valid threshold");
1213        assert_eq!(summary.count, 0);
1214        assert_eq!(summary.mean_energy_mev, 0.0);
1215        assert_eq!(summary.runaway_fraction, 0.0);
1216    }
1217
1218    #[test]
1219    fn test_population_summary_rejects_invalid_threshold() {
1220        let particles = seed_alpha_test_particles(4, 6.2, 0.0, 3.5, 0.2, 1.0).expect("valid seeds");
1221        for bad in [f64::NAN, -0.01] {
1222            let err = summarize_particle_population(&particles, bad)
1223                .expect_err("invalid threshold must error");
1224            match err {
1225                FusionError::PhysicsViolation(msg) => {
1226                    assert!(msg.contains("runaway_threshold_mev"));
1227                }
1228                other => panic!("Unexpected error: {other:?}"),
1229            }
1230        }
1231    }
1232
1233    #[test]
1234    fn test_blend_particle_current_rejects_non_finite_grid_mesh() {
1235        let mut grid = Grid2D::new(8, 8, 1.0, 5.0, -2.0, 2.0);
1236        grid.rr[[0, 0]] = f64::NAN;
1237        let fluid = Array2::from_elem((8, 8), 2.0);
1238        let particle = Array2::from_elem((8, 8), 6.0);
1239        let err = blend_particle_current(&fluid, &particle, &grid, 1.0, 0.5)
1240            .expect_err("non-finite grid mesh must fail");
1241        match err {
1242            FusionError::PhysicsViolation(msg) => {
1243                assert!(msg.contains("mesh coordinates"));
1244            }
1245            other => panic!("Unexpected error: {other:?}"),
1246        }
1247    }
1248
1249    #[test]
1250    fn test_population_summary_rejects_non_finite_particle_state() {
1251        let mut particles =
1252            seed_alpha_test_particles(4, 6.2, 0.0, 3.5, 0.2, 1.0).expect("valid seeds");
1253        particles[0].vx_m_s = f64::INFINITY;
1254        let err = summarize_particle_population(&particles, 1.0)
1255            .expect_err("non-finite particle state must fail");
1256        match err {
1257            FusionError::PhysicsViolation(msg) => {
1258                assert!(msg.contains("velocity"));
1259            }
1260            other => panic!("Unexpected error: {other:?}"),
1261        }
1262    }
1263
1264    #[test]
1265    fn test_toroidal_current_deposition_rejects_invalid_particle_state() {
1266        let grid = Grid2D::new(17, 17, 1.0, 9.0, -4.0, 4.0);
1267        let mut particles =
1268            seed_alpha_test_particles(4, 6.2, 0.0, 3.5, 0.2, 1.0).expect("valid seeds");
1269        particles[1].weight = f64::NAN;
1270        let err = deposit_toroidal_current_density(&particles, &grid)
1271            .expect_err("invalid particle state must fail");
1272        match err {
1273            FusionError::PhysicsViolation(msg) => {
1274                assert!(msg.contains("weight"));
1275            }
1276            other => panic!("Unexpected error: {other:?}"),
1277        }
1278    }
1279
1280    #[test]
1281    fn test_toroidal_current_deposition_rejects_empty_population() {
1282        let grid = Grid2D::new(17, 17, 1.0, 9.0, -4.0, 4.0);
1283        let err = deposit_toroidal_current_density(&[], &grid)
1284            .expect_err("empty particle population must fail");
1285        match err {
1286            FusionError::PhysicsViolation(msg) => {
1287                assert!(msg.contains("non-empty"));
1288            }
1289            other => panic!("Unexpected error: {other:?}"),
1290        }
1291    }
1292
1293    #[test]
1294    fn test_toroidal_current_deposition_supports_descending_axes() {
1295        let grid = Grid2D::new(17, 17, 9.0, 1.0, 4.0, -4.0);
1296        let particles = vec![ChargedParticle {
1297            x_m: 5.0,
1298            y_m: 0.0,
1299            z_m: 0.0,
1300            vx_m_s: 0.0,
1301            vy_m_s: 100_000.0,
1302            vz_m_s: 0.0,
1303            charge_c: 1.602_176_634e-19,
1304            mass_kg: 1.672_621_923_69e-27,
1305            weight: 2.0e16,
1306        }];
1307        let j = deposit_toroidal_current_density(&particles, &grid)
1308            .expect("descending axes should still deposit current");
1309        let sum_abs = j.iter().map(|v| v.abs()).sum::<f64>();
1310        assert!(
1311            sum_abs > 0.0,
1312            "Expected non-zero toroidal deposition on descending axes"
1313        );
1314    }
1315
1316    #[test]
1317    fn test_boris_push_rejects_invalid_runtime_inputs() {
1318        let mut particle =
1319            seed_alpha_test_particles(1, 6.2, 0.0, 3.5, 0.2, 1.0).expect("valid seeds")[0];
1320        let err = boris_push_step(&mut particle, [0.0, 0.0, f64::NAN], [0.0, 0.0, 2.5], 1e-9)
1321            .expect_err("non-finite electric field must fail");
1322        match err {
1323            FusionError::PhysicsViolation(msg) => {
1324                assert!(msg.contains("electric_v_m"));
1325            }
1326            other => panic!("Unexpected error: {other:?}"),
1327        }
1328        let err = boris_push_step(&mut particle, [0.0, 0.0, 0.0], [0.0, 0.0, 2.5], 0.0)
1329            .expect_err("non-positive dt must fail");
1330        match err {
1331            FusionError::PhysicsViolation(msg) => {
1332                assert!(msg.contains("dt_s"));
1333            }
1334            other => panic!("Unexpected error: {other:?}"),
1335        }
1336    }
1337
1338    #[test]
1339    fn test_advance_particles_boris_rejects_invalid_particle_state() {
1340        let mut particles =
1341            seed_alpha_test_particles(2, 6.2, 0.0, 3.5, 0.2, 1.0).expect("valid seeds");
1342        particles[0].mass_kg = 0.0;
1343        let err =
1344            advance_particles_boris(&mut particles, [0.0, 0.0, 0.0], [0.0, 0.0, 2.5], 1e-9, 1)
1345                .expect_err("invalid particle mass must fail");
1346        match err {
1347            FusionError::PhysicsViolation(msg) => {
1348                assert!(msg.contains("particle[0]"));
1349            }
1350            other => panic!("Unexpected error: {other:?}"),
1351        }
1352    }
1353
1354    #[test]
1355    fn test_advance_particles_boris_rejects_empty_slice_and_zero_steps() {
1356        let mut empty: Vec<ChargedParticle> = Vec::new();
1357        let err = advance_particles_boris(&mut empty, [0.0, 0.0, 0.0], [0.0, 0.0, 2.5], 1e-9, 1)
1358            .expect_err("empty particle slice must fail");
1359        match err {
1360            FusionError::PhysicsViolation(msg) => {
1361                assert!(msg.contains("non-empty"));
1362            }
1363            other => panic!("Unexpected error: {other:?}"),
1364        }
1365
1366        let mut particles =
1367            seed_alpha_test_particles(1, 6.2, 0.0, 3.5, 0.2, 1.0).expect("valid seeds");
1368        let err =
1369            advance_particles_boris(&mut particles, [0.0, 0.0, 0.0], [0.0, 0.0, 2.5], 1e-9, 0)
1370                .expect_err("zero steps must fail");
1371        match err {
1372            FusionError::PhysicsViolation(msg) => {
1373                assert!(msg.contains("steps"));
1374            }
1375            other => panic!("Unexpected error: {other:?}"),
1376        }
1377    }
1378
1379    // ═══════════════════════════════════════════════════════════════════
1380    // Coulomb collision tests
1381    // ═══════════════════════════════════════════════════════════════════
1382
1383    fn default_collision_params() -> CoulombCollisionParams {
1384        CoulombCollisionParams {
1385            n_e: 1.0e20,
1386            t_e_kev: 10.0,
1387            t_i_kev: 8.0,
1388            a_i: 2.0, // deuterium
1389            z_i: 1.0,
1390            z_eff: 1.5,
1391        }
1392    }
1393
1394    #[test]
1395    fn test_coulomb_logarithm_range() {
1396        let ln_lam = coulomb_logarithm(1e20, 10.0).unwrap();
1397        assert!((5.0..=30.0).contains(&ln_lam), "ln_lambda={ln_lam}");
1398        // For fusion plasma, NRL formula gives ln Λ ≈ 10-20 range
1399        assert!(
1400            ln_lam > 5.0 && ln_lam < 25.0,
1401            "fusion plasma ln_lambda={ln_lam}"
1402        );
1403    }
1404
1405    #[test]
1406    fn test_coulomb_logarithm_clamp() {
1407        // Very cold dense plasma -> clamped to 5
1408        let low = coulomb_logarithm(1e30, 0.001).unwrap();
1409        assert!((low - 5.0).abs() < 1e-12, "should clamp to 5, got {low}");
1410        // Moderately hot tenuous plasma -> still within [5,30]
1411        let mid = coulomb_logarithm(1e10, 1000.0).unwrap();
1412        assert!(
1413            (5.0..=30.0).contains(&mid),
1414            "should be in [5,30], got {mid}"
1415        );
1416    }
1417
1418    #[test]
1419    fn test_coulomb_logarithm_rejects_invalid() {
1420        assert!(coulomb_logarithm(-1.0, 10.0).is_err());
1421        assert!(coulomb_logarithm(1e20, 0.0).is_err());
1422        assert!(coulomb_logarithm(1e20, f64::NAN).is_err());
1423    }
1424
1425    #[test]
1426    fn test_spitzer_slowing_down_time() {
1427        let ln_lam = coulomb_logarithm(1e20, 10.0).unwrap();
1428        let tau = spitzer_slowing_down_time(ALPHA_MASS_KG, 2.0, 1e20, 10.0, ln_lam).unwrap();
1429        // For 3.5 MeV alpha in ITER plasma, τ_s ~ 0.1-1.0 s
1430        assert!(
1431            tau > 0.01 && tau < 10.0,
1432            "spitzer time {tau} s out of range"
1433        );
1434    }
1435
1436    #[test]
1437    fn test_critical_velocity() {
1438        let v_c = critical_velocity(10.0, 2.0, 1.0, 1.5).unwrap();
1439        // v_c should be much less than speed of light but significant
1440        assert!(v_c > 1e5 && v_c < 1e8, "v_c={v_c} m/s");
1441    }
1442
1443    #[test]
1444    fn test_collision_frequencies_positive() {
1445        let params = default_collision_params();
1446        let ln_lam = coulomb_logarithm(params.n_e, params.t_e_kev).unwrap();
1447        let tau_s =
1448            spitzer_slowing_down_time(ALPHA_MASS_KG, 2.0, params.n_e, params.t_e_kev, ln_lam)
1449                .unwrap();
1450        let v_c = critical_velocity(params.t_e_kev, params.a_i, params.z_i, params.z_eff).unwrap();
1451        let speed = 1e7; // fast alpha
1452        let (nu_s, nu_d, nu_e) = collision_frequencies(speed, &params, ln_lam, tau_s, v_c).unwrap();
1453        assert!(nu_s > 0.0, "slowing-down frequency must be > 0");
1454        assert!(nu_d > 0.0, "deflection frequency must be > 0");
1455        assert!(nu_e > 0.0, "energy diffusion frequency must be > 0");
1456    }
1457
1458    #[test]
1459    fn test_collision_step_slows_alpha() {
1460        let params = default_collision_params();
1461        let mut particles = seed_alpha_test_particles(1, 6.2, 0.0, 3.5, 0.6, 1e12).unwrap();
1462        let e0 = particles[0].kinetic_energy_mev();
1463        let mut rng: u64 = 12345;
1464        // Many small collision steps should reduce energy
1465        for _ in 0..100 {
1466            collision_step(&mut particles[0], &params, 1e-3, &mut rng).unwrap();
1467        }
1468        let e1 = particles[0].kinetic_energy_mev();
1469        assert!(e1 < e0, "collisions should slow alpha: {e0} -> {e1} MeV");
1470    }
1471
1472    #[test]
1473    fn test_batch_collisions_deterministic() {
1474        let params = default_collision_params();
1475        let mut p1 = seed_alpha_test_particles(8, 6.2, 0.0, 3.5, 0.6, 1e12).unwrap();
1476        let mut p2 = p1.clone();
1477        apply_coulomb_collisions(&mut p1, &params, 1e-4, 42).unwrap();
1478        apply_coulomb_collisions(&mut p2, &params, 1e-4, 42).unwrap();
1479        for (a, b) in p1.iter().zip(p2.iter()) {
1480            assert!((a.vx_m_s - b.vx_m_s).abs() < 1e-12, "determinism broken");
1481            assert!((a.vy_m_s - b.vy_m_s).abs() < 1e-12, "determinism broken");
1482            assert!((a.vz_m_s - b.vz_m_s).abs() < 1e-12, "determinism broken");
1483        }
1484    }
1485
1486    #[test]
1487    fn test_collision_rejects_invalid_params() {
1488        let mut bad = default_collision_params();
1489        bad.n_e = -1.0;
1490        let mut p = seed_alpha_test_particles(1, 6.2, 0.0, 3.5, 0.6, 1e12).unwrap();
1491        assert!(apply_coulomb_collisions(&mut p, &bad, 1e-4, 42).is_err());
1492    }
1493
1494    #[test]
1495    fn test_perpendicular_basis_orthonormal() {
1496        for v_hat in [
1497            [1.0, 0.0, 0.0],
1498            [0.0, 1.0, 0.0],
1499            [0.0, 0.0, 1.0],
1500            [0.577, 0.577, 0.577],
1501        ] {
1502            let norm = dot(v_hat, v_hat).sqrt();
1503            let v_hat = [v_hat[0] / norm, v_hat[1] / norm, v_hat[2] / norm];
1504            let (e1, e2) = perpendicular_basis(v_hat);
1505            // e1 ⊥ v_hat
1506            assert!(dot(e1, v_hat).abs() < 1e-10, "e1 not perp to v_hat");
1507            // e2 ⊥ v_hat
1508            assert!(dot(e2, v_hat).abs() < 1e-10, "e2 not perp to v_hat");
1509            // e1 ⊥ e2
1510            assert!(dot(e1, e2).abs() < 1e-10, "e1 not perp to e2");
1511            // |e1| ≈ 1, |e2| ≈ 1
1512            assert!((dot(e1, e1).sqrt() - 1.0).abs() < 1e-10);
1513            assert!((dot(e2, e2).sqrt() - 1.0).abs() < 1e-10);
1514        }
1515    }
1516}