1use crate::pedestal::{PedestalConfig, PedestalModel};
7use fusion_math::tridiag::thomas_solve;
8use fusion_types::error::{FusionError, FusionResult};
9use fusion_types::state::RadialProfiles;
10use ndarray::Array1;
11
12const TRANSPORT_NR: usize = 50;
14
15const CHI_BASE_DEFAULT: f64 = 0.5;
17
18const CHI_NC_FLOOR: f64 = 0.01;
20
21const NC_BOLTZMANN_J_PER_KEV: f64 = 1.602_176_634e-16;
23const NC_ELEM_CHARGE: f64 = 1.602_176_634e-19;
25const NC_PROTON_MASS: f64 = 1.672_621_923_69e-27;
27
28const CHI_TURB: f64 = 5.0;
30
31const D_IMPURITY: f64 = 1.0;
33
34const CRIT_GRADIENT: f64 = 2.0;
36
37const HMODE_CHI_MIN_FACTOR: f64 = 0.08;
39
40const HMODE_POWER_THRESHOLD: f64 = 30.0;
42
43const EDGE_TEMPERATURE: f64 = 0.1;
45const MAX_TEMPERATURE: f64 = 100.0;
47
48const COOLING_FACTOR: f64 = 5.0;
50
51const HEATING_WIDTH: f64 = 0.1;
53
54const TOROIDAL_COUPLING_MAX_FACTOR: f64 = 3.0;
56
57#[derive(Debug, Clone, PartialEq)]
59pub struct NeoclassicalParams {
60 pub r_major: f64,
62 pub a_minor: f64,
64 pub b_toroidal: f64,
66 pub a_ion: f64,
68 pub z_eff: f64,
70 pub q_profile: Array1<f64>,
72}
73
74impl Default for NeoclassicalParams {
75 fn default() -> Self {
76 Self {
77 r_major: 6.2,
78 a_minor: 2.0,
79 b_toroidal: 5.3,
80 a_ion: 2.0,
81 z_eff: 1.5,
82 q_profile: Array1::linspace(1.0, 3.0, TRANSPORT_NR),
83 }
84 }
85}
86
87pub fn chang_hinton_chi(
94 rho: f64,
95 t_i_kev: f64,
96 n_e_19: f64,
97 q: f64,
98 params: &NeoclassicalParams,
99) -> f64 {
100 if rho <= 0.0 || rho > 1.0 || t_i_kev <= 0.0 || n_e_19 <= 0.0 || q <= 0.0 {
101 return CHI_NC_FLOOR;
102 }
103 let epsilon = (rho * params.a_minor) / params.r_major;
104 if epsilon <= 1e-6 {
105 return CHI_NC_FLOOR;
106 }
107
108 let t_i_j = t_i_kev * NC_BOLTZMANN_J_PER_KEV;
109 let m_i = params.a_ion * NC_PROTON_MASS;
110 let n_e = n_e_19 * 1e19;
111
112 let v_ti = (2.0 * t_i_j / m_i).sqrt();
114 if !v_ti.is_finite() || v_ti <= 0.0 {
115 return CHI_NC_FLOOR;
116 }
117
118 let rho_i = m_i * v_ti / (NC_ELEM_CHARGE * params.b_toroidal);
120
121 let ln_lambda = 17.0; let eps0 = 8.854_187_812_8e-12;
125 let e4 = NC_ELEM_CHARGE.powi(4);
126 let nu_ii = n_e * params.z_eff * params.z_eff * e4 * ln_lambda
127 / (12.0 * std::f64::consts::PI.powf(1.5) * eps0 * eps0 * m_i.sqrt() * t_i_j.powf(1.5));
128
129 let eps32 = epsilon.powf(1.5);
131 let nu_star = nu_ii * q * params.r_major / (eps32 * v_ti);
132
133 let alpha = epsilon; let chi = 0.66 * (1.0 + 1.54 * alpha) * q * q * rho_i * rho_i * nu_ii
138 / (eps32 * (1.0 + 0.74 * nu_star.powf(2.0 / 3.0)));
139
140 if chi.is_finite() && chi > 0.0 {
141 chi.max(CHI_NC_FLOOR)
142 } else {
143 CHI_NC_FLOOR
144 }
145}
146
147pub fn chang_hinton_chi_profile(
149 rho: &Array1<f64>,
150 t_i_kev: &Array1<f64>,
151 n_e_19: &Array1<f64>,
152 q: &Array1<f64>,
153 params: &NeoclassicalParams,
154) -> Array1<f64> {
155 if rho.len() != t_i_kev.len() || rho.len() != n_e_19.len() || rho.len() != q.len() {
156 return Array1::from_elem(rho.len(), CHI_NC_FLOOR);
157 }
158 Array1::from_iter(
159 (0..rho.len()).map(|i| chang_hinton_chi(rho[i], t_i_kev[i], n_e_19[i], q[i], params)),
160 )
161}
162
163pub fn gyro_bohm_chi_profile(
165 rho: &Array1<f64>,
166 t_i_kev: &Array1<f64>,
167 b_field: f64,
168 mass_amu: f64,
169 chi_gb_coeff: f64,
170) -> Array1<f64> {
171 let n = rho.len();
172 if n == 0 || t_i_kev.len() != n {
173 return Array1::zeros(n);
174 }
175 if !b_field.is_finite() || b_field.abs() < 1e-12 || !mass_amu.is_finite() || mass_amu <= 0.0 {
176 return Array1::from_elem(n, CHI_NC_FLOOR);
177 }
178 let b2 = b_field * b_field;
179 let mass_term = mass_amu.sqrt();
180 let coeff = chi_gb_coeff.abs().max(1e-9);
181 Array1::from_iter((0..n).map(|i| {
182 let t = t_i_kev[i].max(0.01);
183 let r = rho[i].clamp(0.0, 1.0);
184 let chi = coeff * (t.powf(1.5) / b2) * (1.0 + r) / mass_term;
185 if chi.is_finite() && chi > 0.0 {
186 chi.max(CHI_NC_FLOOR)
187 } else {
188 CHI_NC_FLOOR
189 }
190 }))
191}
192
193pub fn sauter_bootstrap_current_profile(
195 rho: &Array1<f64>,
196 t_e_kev: &Array1<f64>,
197 t_i_kev: &Array1<f64>,
198 n_e_19: &Array1<f64>,
199 q: &Array1<f64>,
200 epsilon: &Array1<f64>,
201 b_field: f64,
202) -> Array1<f64> {
203 let n = rho.len();
204 if n < 3
205 || t_e_kev.len() != n
206 || t_i_kev.len() != n
207 || n_e_19.len() != n
208 || q.len() != n
209 || epsilon.len() != n
210 || !b_field.is_finite()
211 || b_field <= 0.0
212 {
213 return Array1::zeros(n);
214 }
215
216 let e_charge = NC_ELEM_CHARGE;
217 let m_e = 9.109_383_70e-31;
218 let eps0: f64 = 8.854_187_812e-12;
219 let z_eff = 1.5;
220 let mut j_bs = Array1::zeros(n);
221
222 for i in 1..(n - 1) {
223 let eps = epsilon[i];
224 let te = t_e_kev[i];
225 let ti = t_i_kev[i];
226 let ne19 = n_e_19[i];
227 let qi = q[i];
228 if eps < 1e-6 || te <= 0.0 || ti <= 0.0 || ne19 <= 0.0 || qi <= 0.0 {
229 continue;
230 }
231
232 let eps_sq = (eps * eps).clamp(0.0, 0.999_999);
233 let trapped_denom = (1.0 - eps_sq).sqrt() * (1.0 + 1.46 * eps.sqrt());
234 if trapped_denom <= 1e-12 {
235 continue;
236 }
237 let mut f_t = 1.0 - (1.0 - eps).powi(2) / trapped_denom;
238 f_t = f_t.clamp(0.0, 1.0);
239
240 let t_e_j = te * 1e3 * e_charge;
241 if !t_e_j.is_finite() || t_e_j <= 0.0 {
242 continue;
243 }
244 let v_te = (2.0 * t_e_j / m_e).sqrt();
245 if !v_te.is_finite() || v_te <= 0.0 {
246 continue;
247 }
248
249 let n_e = ne19 * 1e19;
250 let ln_lambda = 17.0;
251 let nu_ei = n_e * z_eff * e_charge.powi(4) * ln_lambda
252 / (12.0 * std::f64::consts::PI.powf(1.5) * eps0.powi(2) * m_e.sqrt() * t_e_j.powf(1.5));
253 if !nu_ei.is_finite() || nu_ei <= 0.0 {
254 continue;
255 }
256 let nu_star_e = nu_ei * qi / (eps.powf(1.5).max(1e-9) * v_te);
257 let alpha_31 = 1.0 / (1.0 + 0.36 / z_eff);
258 let l31 = f_t * alpha_31
259 / (1.0 + alpha_31 * nu_star_e.sqrt() + 0.25 * nu_star_e * (1.0 - f_t).powi(2));
260 let mut l32 = f_t * (0.05 + 0.62 * z_eff) / (z_eff * (1.0 + 0.44 * z_eff));
261 l32 /= 1.0 + 0.22 * nu_star_e.sqrt() + 0.19 * nu_star_e * (1.0 - f_t);
262 let l34 = l31 * ti / te.max(0.01);
263
264 let dr = (rho[i + 1] - rho[i - 1]).abs().max(1e-12);
265 let dn_dr = (n_e_19[i + 1] - n_e_19[i - 1]) * 1e19 / dr;
266 let dte_dr = (t_e_kev[i + 1] - t_e_kev[i - 1]) * 1e3 * e_charge / dr;
267 let dti_dr = (t_i_kev[i + 1] - t_i_kev[i - 1]) * 1e3 * e_charge / dr;
268
269 let b_pol = (b_field * eps / qi.max(0.1)).max(1e-10);
270 let p_e = n_e * t_e_j;
271 let denom_ti = (ti * 1e3 * e_charge).max(1e-30);
272 let j_val = -(p_e / b_pol)
273 * (l31 * dn_dr / n_e.max(1e10)
274 + l32 * dte_dr / t_e_j.max(1e-30)
275 + l34 * dti_dr / denom_ti);
276 if j_val.is_finite() {
277 j_bs[i] = j_val;
278 }
279 }
280 j_bs
281}
282
283pub fn build_cn_tridiag(
289 chi: &Array1<f64>,
290 rho: &Array1<f64>,
291 dt: f64,
292) -> (Vec<f64>, Vec<f64>, Vec<f64>, Vec<f64>) {
293 let n = rho.len();
294 let mut a = vec![0.0; n];
295 let mut b = vec![1.0; n];
296 let mut c = vec![0.0; n];
297 let mut d = vec![0.0; n];
298 if n < 3 || chi.len() != n || !dt.is_finite() || dt <= 0.0 {
299 return (a, b, c, d);
300 }
301
302 let dr = (rho[1] - rho[0]).abs().max(1e-12);
303 for i in 1..(n - 1) {
304 let r = rho[i].abs().max(1e-6);
305 let chi_i = chi[i].max(0.0);
306 let chi_ip = 0.5 * (chi_i + chi[i + 1].max(0.0));
307 let chi_im = 0.5 * (chi_i + chi[i - 1].max(0.0));
308 let r_ip = r + 0.5 * dr;
309 let r_im = (r - 0.5 * dr).max(1e-6);
310 let coeff_ip = chi_ip * r_ip / (r * dr * dr);
311 let coeff_im = chi_im * r_im / (r * dr * dr);
312
313 b[i] = 1.0 + 0.5 * dt * (coeff_ip + coeff_im);
314 c[i] = -0.5 * dt * coeff_ip;
315 a[i] = -0.5 * dt * coeff_im;
316 }
317 d[0] = EDGE_TEMPERATURE;
318 d[n - 1] = EDGE_TEMPERATURE;
319 (a, b, c, d)
320}
321
322pub fn transport_step(solver: &mut TransportSolver, p_aux_mw: f64, dt: f64) -> FusionResult<()> {
324 if !p_aux_mw.is_finite() || p_aux_mw < 0.0 {
325 return Err(FusionError::ConfigError(format!(
326 "transport step requires finite p_aux_mw >= 0, got {p_aux_mw}"
327 )));
328 }
329 if !dt.is_finite() || dt <= 0.0 {
330 return Err(FusionError::ConfigError(format!(
331 "transport step requires finite dt > 0, got {dt}"
332 )));
333 }
334
335 let dt_prev = solver.dt;
336 solver.dt = dt;
337 let step_result = (|| -> FusionResult<()> {
338 solver.update_transport_model(p_aux_mw)?;
339 let te_before = solver.profiles.te.clone();
340 solver.evolve_profiles(p_aux_mw)?;
341
342 let (a, b, c, mut d) = build_cn_tridiag(&solver.chi, &solver.profiles.rho, dt);
344 let n = solver.profiles.te.len();
345 if n >= 3 {
346 for (i, d_i) in d.iter_mut().enumerate().skip(1).take(n - 2) {
347 *d_i = solver.profiles.te[i];
348 }
349 let solved = thomas_solve(&a, &b, &c, &d);
350 for i in 1..(n - 1) {
351 let val = solved[i];
352 solver.profiles.te[i] = if val.is_finite() {
353 val.clamp(EDGE_TEMPERATURE, MAX_TEMPERATURE)
354 } else {
355 te_before[i]
356 };
357 solver.profiles.ti[i] = solver.profiles.te[i];
358 }
359 solver.profiles.te[n - 1] = EDGE_TEMPERATURE;
360 solver.profiles.ti[n - 1] = EDGE_TEMPERATURE;
361 }
362
363 if solver.profiles.te.iter().any(|v| !v.is_finite())
364 || solver.profiles.ti.iter().any(|v| !v.is_finite())
365 || solver.profiles.ne.iter().any(|v| !v.is_finite())
366 {
367 return Err(FusionError::ConfigError(
368 "transport step produced non-finite profile outputs".to_string(),
369 ));
370 }
371 Ok(())
372 })();
373 solver.dt = dt_prev;
374 step_result
375}
376
377pub struct TransportSolver {
379 pub profiles: RadialProfiles,
380 pub chi: Array1<f64>,
381 pub dt: f64,
382 pub pedestal: PedestalModel,
383 pub toroidal_mode_amplitudes: Vec<f64>,
384 pub toroidal_coupling_gain: f64,
385 pub neoclassical: Option<NeoclassicalParams>,
387}
388
389impl TransportSolver {
390 pub fn new() -> Self {
392 let rho = Array1::linspace(0.0, 1.0, TRANSPORT_NR);
393
394 let te = Array1::from_shape_fn(TRANSPORT_NR, |i| {
396 let r: f64 = rho[i];
397 10.0 * (1.0 - r * r).max(0.0) + EDGE_TEMPERATURE
398 });
399 let ti = te.clone();
400 let ne = Array1::from_shape_fn(TRANSPORT_NR, |i| {
401 let r: f64 = rho[i];
402 10.0 * (1.0 - r * r).max(0.0) + 0.5
403 });
404 let n_impurity = Array1::zeros(TRANSPORT_NR);
405 let chi = Array1::from_elem(TRANSPORT_NR, CHI_BASE_DEFAULT);
406 let pedestal = PedestalModel::new(PedestalConfig {
407 beta_p_ped: 0.35,
408 rho_s: 2.0e-3,
409 r_major: 6.2,
410 alpha_crit: 2.5,
411 tau_elm: 1.0e-3,
412 })
413 .expect("default pedestal config must be valid");
414
415 TransportSolver {
416 profiles: RadialProfiles {
417 rho,
418 te,
419 ti,
420 ne,
421 n_impurity,
422 },
423 chi,
424 dt: 0.01,
425 pedestal,
426 toroidal_mode_amplitudes: vec![0.0; 3],
427 toroidal_coupling_gain: 0.0,
428 neoclassical: None,
429 }
430 }
431
432 pub fn set_neoclassical(&mut self, params: NeoclassicalParams) -> FusionResult<()> {
436 if !params.r_major.is_finite() || params.r_major <= 0.0 {
437 return Err(FusionError::PhysicsViolation(
438 "neoclassical r_major must be finite and > 0".into(),
439 ));
440 }
441 if !params.a_minor.is_finite() || params.a_minor <= 0.0 {
442 return Err(FusionError::PhysicsViolation(
443 "neoclassical a_minor must be finite and > 0".into(),
444 ));
445 }
446 if !params.b_toroidal.is_finite() || params.b_toroidal <= 0.0 {
447 return Err(FusionError::PhysicsViolation(
448 "neoclassical b_toroidal must be finite and > 0".into(),
449 ));
450 }
451 if !params.a_ion.is_finite() || params.a_ion <= 0.0 {
452 return Err(FusionError::PhysicsViolation(
453 "neoclassical a_ion must be finite and > 0".into(),
454 ));
455 }
456 if !params.z_eff.is_finite() || params.z_eff < 1.0 {
457 return Err(FusionError::PhysicsViolation(
458 "neoclassical z_eff must be finite and >= 1".into(),
459 ));
460 }
461 if params.q_profile.len() != self.profiles.rho.len() {
462 return Err(FusionError::ConfigError(format!(
463 "neoclassical q_profile length {} must match radial grid {}",
464 params.q_profile.len(),
465 self.profiles.rho.len()
466 )));
467 }
468 if params.q_profile.iter().any(|v| !v.is_finite() || *v <= 0.0) {
469 return Err(FusionError::PhysicsViolation(
470 "neoclassical q_profile values must be finite and > 0".into(),
471 ));
472 }
473 self.neoclassical = Some(params);
474 Ok(())
475 }
476
477 pub fn set_toroidal_mode_spectrum(
482 &mut self,
483 amplitudes: &[f64],
484 gain: f64,
485 ) -> FusionResult<()> {
486 if !gain.is_finite() || gain < 0.0 {
487 return Err(FusionError::PhysicsViolation(
488 "toroidal coupling gain must be finite and >= 0".to_string(),
489 ));
490 }
491 if amplitudes.iter().any(|a| !a.is_finite() || *a < 0.0) {
492 return Err(FusionError::PhysicsViolation(
493 "toroidal mode amplitudes must be finite and >= 0".to_string(),
494 ));
495 }
496 self.toroidal_mode_amplitudes.clear();
497 self.toroidal_mode_amplitudes
498 .extend(amplitudes.iter().copied());
499 self.toroidal_coupling_gain = gain;
500 Ok(())
501 }
502
503 fn toroidal_mode_coupling_factor(&self, rho: f64) -> f64 {
504 if self.toroidal_coupling_gain <= 0.0 || self.toroidal_mode_amplitudes.is_empty() {
505 return 1.0;
506 }
507
508 let spectral_rms = self
510 .toroidal_mode_amplitudes
511 .iter()
512 .enumerate()
513 .map(|(idx, amp)| {
514 let n = (idx + 1) as f64;
515 (n * amp).powi(2)
516 })
517 .sum::<f64>()
518 .sqrt();
519 if spectral_rms <= 1e-12 {
520 return 1.0;
521 }
522
523 let envelope = rho.clamp(0.0, 1.0).powi(2);
525 (1.0 + self.toroidal_coupling_gain * envelope * spectral_rms)
526 .clamp(1.0, TOROIDAL_COUPLING_MAX_FACTOR)
527 }
528
529 pub fn update_transport_model(&mut self, p_aux_mw: f64) -> FusionResult<()> {
534 if !p_aux_mw.is_finite() || p_aux_mw < 0.0 {
535 return Err(FusionError::ConfigError(format!(
536 "transport update requires finite p_aux_mw >= 0, got {p_aux_mw}"
537 )));
538 }
539 let n = self.profiles.rho.len();
540 if n < 2 {
541 return Err(FusionError::ConfigError(
542 "transport update requires at least 2 radial points".to_string(),
543 ));
544 }
545 let dr = if n > 1 { 1.0 / (n as f64 - 1.0) } else { 1.0 };
546 if !dr.is_finite() || dr <= 0.0 {
547 return Err(FusionError::ConfigError(format!(
548 "transport update produced invalid dr={dr}"
549 )));
550 }
551 let ped_width = self.pedestal.pedestal_width();
552 if !ped_width.is_finite() || ped_width <= 0.0 {
553 return Err(FusionError::ConfigError(format!(
554 "transport pedestal width must be finite and > 0, got {ped_width}"
555 )));
556 }
557 let barrier_rho = (1.0 - ped_width).clamp(0.75, 0.995);
558
559 for i in 0..n {
560 if !self.profiles.te[i].is_finite() || !self.profiles.rho[i].is_finite() {
561 return Err(FusionError::ConfigError(format!(
562 "transport profile contains non-finite values at index {}",
563 i
564 )));
565 }
566 let grad_t = if i == 0 {
568 (self.profiles.te[1] - self.profiles.te[0]) / dr
569 } else if i == n - 1 {
570 (self.profiles.te[n - 1] - self.profiles.te[n - 2]) / dr
571 } else {
572 (self.profiles.te[i + 1] - self.profiles.te[i - 1]) / (2.0 * dr)
573 };
574 if !grad_t.is_finite() {
575 return Err(FusionError::ConfigError(format!(
576 "transport gradient became non-finite at index {}",
577 i
578 )));
579 }
580
581 let excess = (-grad_t - CRIT_GRADIENT).max(0.0);
583 let chi_base = if let Some(ref nc) = self.neoclassical {
584 let rho_val = self.profiles.rho[i];
585 let t_i = self.profiles.ti[i].max(0.01);
586 let n_e_19 = self.profiles.ne[i].max(0.01);
587 let q = nc.q_profile[i];
588 chang_hinton_chi(rho_val, t_i, n_e_19, q, nc)
589 } else {
590 CHI_BASE_DEFAULT
591 };
592 self.chi[i] = chi_base + CHI_TURB * excess;
593
594 self.chi[i] *= self.toroidal_mode_coupling_factor(self.profiles.rho[i]);
596
597 if p_aux_mw > HMODE_POWER_THRESHOLD && self.profiles.rho[i] >= barrier_rho {
599 let edge_weight =
600 ((self.profiles.rho[i] - barrier_rho) / ped_width.max(1e-5)).clamp(0.0, 1.0);
601 let factor = (1.0 - 0.92 * edge_weight).clamp(HMODE_CHI_MIN_FACTOR, 1.0);
602 self.chi[i] *= factor;
603 }
604 if !self.chi[i].is_finite() {
605 return Err(FusionError::ConfigError(format!(
606 "transport chi became non-finite at index {}",
607 i
608 )));
609 }
610 }
611 Ok(())
612 }
613
614 fn compute_pedestal_pressure_gradient(&self) -> FusionResult<f64> {
615 let n = self.profiles.rho.len();
616 if n < 3 {
617 return Err(FusionError::ConfigError(
618 "transport pressure-gradient calculation requires at least 3 radial points"
619 .to_string(),
620 ));
621 }
622
623 let dr = 1.0 / (n as f64 - 1.0);
624 let ped_width = self.pedestal.pedestal_width();
625 let barrier_rho = (1.0 - ped_width).clamp(0.75, 0.995);
626 let pressure =
627 &self.profiles.ne * (&self.profiles.te + &self.profiles.ti).mapv(|v| v.max(0.0));
628
629 let mut max_grad: f64 = 0.0;
630 for i in 1..n - 1 {
631 if self.profiles.rho[i] < barrier_rho {
632 continue;
633 }
634 let grad = (pressure[i + 1] - pressure[i - 1]) / (2.0 * dr);
635 max_grad = max_grad.max(grad.abs());
636 }
637 if !max_grad.is_finite() {
638 return Err(FusionError::ConfigError(
639 "transport pressure-gradient calculation produced non-finite result".to_string(),
640 ));
641 }
642 Ok(max_grad)
643 }
644
645 pub fn evolve_profiles(&mut self, p_aux_mw: f64) -> FusionResult<()> {
649 if !p_aux_mw.is_finite() || p_aux_mw < 0.0 {
650 return Err(FusionError::ConfigError(format!(
651 "transport evolve requires finite p_aux_mw >= 0, got {p_aux_mw}"
652 )));
653 }
654 let n = self.profiles.rho.len();
655 if n < 2 {
656 return Err(FusionError::ConfigError(
657 "transport evolve requires at least 2 radial points".to_string(),
658 ));
659 }
660 let dr = if n > 1 { 1.0 / (n as f64 - 1.0) } else { 1.0 };
661 let dt = self.dt;
662 if !dt.is_finite() || dt <= 0.0 {
663 return Err(FusionError::ConfigError(format!(
664 "transport evolve requires finite dt > 0, got {dt}"
665 )));
666 }
667
668 let te_old = self.profiles.te.clone();
669
670 for i in 1..n - 1 {
671 let rho_i = self.profiles.rho[i];
672 if !rho_i.is_finite() || rho_i <= 0.0 {
673 return Err(FusionError::ConfigError(format!(
674 "transport evolve requires finite rho > 0 at interior index {}, got {}",
675 i, rho_i
676 )));
677 }
678
679 let flux_plus =
681 0.5 * (self.chi[i] + self.chi[i + 1]) * (te_old[i + 1] - te_old[i]) / dr;
682 let flux_minus =
683 0.5 * (self.chi[i - 1] + self.chi[i]) * (te_old[i] - te_old[i - 1]) / dr;
684 let rho_plus = rho_i + 0.5 * dr;
685 let rho_minus = rho_i - 0.5 * dr;
686 if rho_plus <= 0.0 || rho_minus <= 0.0 {
687 return Err(FusionError::ConfigError(format!(
688 "transport evolve requires positive rho +/- dr/2 at index {}",
689 i
690 )));
691 }
692 let div_flux = (rho_plus * flux_plus - rho_minus * flux_minus) / (rho_i * dr);
693 if !div_flux.is_finite() {
694 return Err(FusionError::ConfigError(format!(
695 "transport evolve produced non-finite div_flux at index {}",
696 i
697 )));
698 }
699
700 let s_heat = p_aux_mw * (-self.profiles.rho[i].powi(2) / HEATING_WIDTH).exp();
702 if !s_heat.is_finite() {
703 return Err(FusionError::ConfigError(format!(
704 "transport evolve produced non-finite heating source at index {}",
705 i
706 )));
707 }
708
709 let s_rad = COOLING_FACTOR
711 * self.profiles.ne[i]
712 * self.profiles.n_impurity[i]
713 * te_old[i].abs().sqrt();
714 if !s_rad.is_finite() {
715 return Err(FusionError::ConfigError(format!(
716 "transport evolve produced non-finite radiation sink at index {}",
717 i
718 )));
719 }
720
721 let te_new = te_old[i] + dt * (div_flux + s_heat - s_rad);
723 if !te_new.is_finite() {
724 return Err(FusionError::ConfigError(format!(
725 "transport evolve produced non-finite temperature at index {}",
726 i
727 )));
728 }
729 self.profiles.te[i] = te_new.clamp(EDGE_TEMPERATURE, MAX_TEMPERATURE);
730 }
731
732 for i in 1..n - 1 {
734 self.profiles.ti[i] = self.profiles.te[i];
735 }
736
737 self.profiles.te[n - 1] = EDGE_TEMPERATURE;
739 self.profiles.ti[n - 1] = EDGE_TEMPERATURE;
740 if self.profiles.te.iter().any(|v| !v.is_finite())
741 || self.profiles.ti.iter().any(|v| !v.is_finite())
742 {
743 return Err(FusionError::ConfigError(
744 "transport evolve produced non-finite profile outputs".to_string(),
745 ));
746 }
747 Ok(())
748 }
749
750 pub fn inject_impurities(&mut self, erosion_rate: f64) -> FusionResult<()> {
752 if !erosion_rate.is_finite() || erosion_rate < 0.0 {
753 return Err(FusionError::ConfigError(format!(
754 "transport impurity injection requires finite erosion_rate >= 0, got {erosion_rate}"
755 )));
756 }
757 if !self.dt.is_finite() || self.dt <= 0.0 {
758 return Err(FusionError::ConfigError(format!(
759 "transport impurity injection requires finite dt > 0, got {}",
760 self.dt
761 )));
762 }
763 let n = self.profiles.rho.len();
764 let dr = if n > 1 { 1.0 / (n as f64 - 1.0) } else { 1.0 };
765
766 if n > 1 {
768 self.profiles.n_impurity[n - 1] += erosion_rate * 1e-18 * self.dt;
769 }
770
771 let n_imp_old = self.profiles.n_impurity.clone();
773 for i in 1..n - 1 {
774 let laplacian = (n_imp_old[i + 1] - 2.0 * n_imp_old[i] + n_imp_old[i - 1]) / (dr * dr);
775 self.profiles.n_impurity[i] =
776 (n_imp_old[i] + D_IMPURITY * self.dt * laplacian).max(0.0);
777 }
778 if self.profiles.n_impurity.iter().any(|v| !v.is_finite()) {
779 return Err(FusionError::ConfigError(
780 "transport impurity injection produced non-finite values".to_string(),
781 ));
782 }
783 Ok(())
784 }
785
786 pub fn step(&mut self, p_aux_mw: f64, erosion_rate: f64) -> FusionResult<()> {
788 if !p_aux_mw.is_finite() || p_aux_mw < 0.0 {
789 return Err(FusionError::ConfigError(format!(
790 "transport step requires finite p_aux_mw >= 0, got {p_aux_mw}"
791 )));
792 }
793 if !erosion_rate.is_finite() || erosion_rate < 0.0 {
794 return Err(FusionError::ConfigError(format!(
795 "transport step requires finite erosion_rate >= 0, got {erosion_rate}"
796 )));
797 }
798 self.pedestal.advance(self.dt);
799 self.update_transport_model(p_aux_mw)?;
800 self.evolve_profiles(p_aux_mw)?;
801
802 if p_aux_mw > HMODE_POWER_THRESHOLD {
803 let grad_p = self.compute_pedestal_pressure_gradient()?;
804 self.pedestal.record_gradient(grad_p);
805 if self.pedestal.is_elm_triggered(grad_p) {
806 self.pedestal.apply_elm_crash(&mut self.profiles);
807 }
808 }
809
810 let n = self.profiles.te.len();
812 if n > 0 {
813 self.profiles.te[n - 1] = EDGE_TEMPERATURE;
814 self.profiles.ti[n - 1] = EDGE_TEMPERATURE;
815 }
816
817 self.inject_impurities(erosion_rate)?;
818 if self.profiles.te.iter().any(|v| !v.is_finite())
819 || self.profiles.ti.iter().any(|v| !v.is_finite())
820 || self.profiles.n_impurity.iter().any(|v| !v.is_finite())
821 || self.chi.iter().any(|v| !v.is_finite())
822 {
823 return Err(FusionError::ConfigError(
824 "transport step produced non-finite runtime state".to_string(),
825 ));
826 }
827 Ok(())
828 }
829}
830
831impl Default for TransportSolver {
832 fn default() -> Self {
833 Self::new()
834 }
835}
836
837#[cfg(test)]
838mod tests {
839 use super::*;
840
841 #[test]
842 fn test_transport_creation() {
843 let ts = TransportSolver::new();
844 assert_eq!(ts.profiles.rho.len(), 50);
845 assert!(ts.profiles.te[0] > ts.profiles.te[49]);
846 }
847
848 #[test]
849 fn test_transport_step_no_panic() {
850 let mut ts = TransportSolver::new();
851 ts.step(50.0, 1e14).expect("valid transport step");
852 assert!(ts.profiles.te.iter().all(|v| v.is_finite()));
854 assert!(ts.profiles.ti.iter().all(|v| v.is_finite()));
855 }
856
857 #[test]
858 fn test_hmode_barrier() {
859 let mut ts = TransportSolver::new();
860
861 ts.update_transport_model(10.0)
863 .expect("valid transport-model update");
864 let chi_edge_no_hmode = ts.chi[48]; ts.update_transport_model(50.0)
868 .expect("valid transport-model update");
869 let chi_edge_hmode = ts.chi[48];
870
871 assert!(
872 chi_edge_hmode < chi_edge_no_hmode,
873 "H-mode should reduce edge chi: {} vs {}",
874 chi_edge_hmode,
875 chi_edge_no_hmode
876 );
877 }
878
879 #[test]
880 fn test_toroidal_mode_coupling_disabled_by_default() {
881 let mut baseline = TransportSolver::new();
882 baseline
883 .update_transport_model(20.0)
884 .expect("valid transport-model update");
885 let baseline_chi = baseline.chi.clone();
886
887 let mut coupled = TransportSolver::new();
888 coupled
889 .set_toroidal_mode_spectrum(&[0.2, 0.1, 0.05], 0.0)
890 .expect("valid spectrum");
891 coupled
892 .update_transport_model(20.0)
893 .expect("valid transport-model update");
894
895 let max_err = baseline_chi
896 .iter()
897 .zip(coupled.chi.iter())
898 .map(|(a, b)| (a - b).abs())
899 .fold(0.0, f64::max);
900 assert!(
901 max_err < 1e-12,
902 "Expected baseline parity when toroidal coupling gain is zero, max_err={max_err}"
903 );
904 }
905
906 #[test]
907 fn test_toroidal_mode_coupling_boosts_edge_more_than_core() {
908 let mut baseline = TransportSolver::new();
909 baseline
910 .update_transport_model(20.0)
911 .expect("valid transport-model update");
912 let base_edge = baseline.chi[48];
913 let base_core = baseline.chi[4];
914
915 let mut coupled = TransportSolver::new();
916 coupled
917 .set_toroidal_mode_spectrum(&[0.2, 0.1, 0.04], 0.7)
918 .expect("valid spectrum");
919 coupled
920 .update_transport_model(20.0)
921 .expect("valid transport-model update");
922 let edge_gain = coupled.chi[48] / base_edge.max(1e-12);
923 let core_gain = coupled.chi[4] / base_core.max(1e-12);
924
925 assert!(
926 edge_gain > 1.0,
927 "Expected edge transport gain > 1 from toroidal coupling, got {edge_gain}"
928 );
929 assert!(
930 edge_gain > core_gain,
931 "Expected stronger edge gain than core gain: edge={edge_gain}, core={core_gain}"
932 );
933 }
934
935 #[test]
936 fn test_toroidal_mode_coupling_factor_is_clamped() {
937 let mut ts = TransportSolver::new();
938 ts.set_toroidal_mode_spectrum(&[10.0, 10.0, 10.0], 10.0)
939 .expect("valid spectrum");
940 let factor = ts.toroidal_mode_coupling_factor(1.0);
941 assert!(
942 (factor - TOROIDAL_COUPLING_MAX_FACTOR).abs() < 1e-12,
943 "Expected clamped factor={}, got {}",
944 TOROIDAL_COUPLING_MAX_FACTOR,
945 factor
946 );
947 }
948
949 #[test]
950 fn test_toroidal_mode_spectrum_rejects_invalid_inputs() {
951 let mut ts = TransportSolver::new();
952 for bad_gain in [f64::NAN, -0.1] {
953 let err = ts
954 .set_toroidal_mode_spectrum(&[0.1, 0.2], bad_gain)
955 .expect_err("invalid gain must error");
956 match err {
957 FusionError::PhysicsViolation(msg) => {
958 assert!(msg.contains("gain"));
959 }
960 other => panic!("Unexpected error: {other:?}"),
961 }
962 }
963 for bad_amp in [f64::NAN, -0.2] {
964 let err = ts
965 .set_toroidal_mode_spectrum(&[0.1, bad_amp], 0.3)
966 .expect_err("invalid amplitude must error");
967 match err {
968 FusionError::PhysicsViolation(msg) => {
969 assert!(msg.contains("amplitudes"));
970 }
971 other => panic!("Unexpected error: {other:?}"),
972 }
973 }
974 }
975
976 #[test]
977 fn test_transport_edge_boundary() {
978 let mut ts = TransportSolver::new();
979 for _ in 0..10 {
980 ts.step(50.0, 0.0).expect("valid transport step");
981 }
982 let last = ts.profiles.te.len() - 1;
984 assert!(
985 (ts.profiles.te[last] - EDGE_TEMPERATURE).abs() < 1e-10,
986 "Edge temperature should be fixed at {EDGE_TEMPERATURE}"
987 );
988 }
989
990 #[test]
991 fn test_elm_crash_reduces_pedestal_temperature() {
992 let mut ts = TransportSolver::new();
993
994 for i in 0..ts.profiles.rho.len() {
996 if ts.profiles.rho[i] < 0.9 {
997 ts.profiles.te[i] = 2.0;
998 ts.profiles.ti[i] = 2.0;
999 ts.profiles.ne[i] = 6.0;
1000 } else {
1001 ts.profiles.te[i] = 9.0;
1002 ts.profiles.ti[i] = 9.0;
1003 ts.profiles.ne[i] = 9.5;
1004 }
1005 }
1006
1007 let edge_idx = ts.profiles.rho.len() - 2;
1008 let te_before = ts.profiles.te[edge_idx];
1009 ts.step(60.0, 0.0).expect("valid transport step");
1010 let te_after = ts.profiles.te[edge_idx];
1011
1012 assert!(
1013 te_after < te_before,
1014 "ELM crash should reduce pedestal Te: before={te_before}, after={te_after}"
1015 );
1016 assert!(ts.pedestal.last_gradient() > 0.0);
1017 }
1018
1019 #[test]
1020 fn test_transport_rejects_invalid_runtime_inputs() {
1021 let mut ts = TransportSolver::new();
1022 let err = ts
1023 .step(f64::NAN, 0.0)
1024 .expect_err("non-finite p_aux must fail");
1025 match err {
1026 FusionError::ConfigError(msg) => assert!(msg.contains("p_aux")),
1027 other => panic!("Unexpected error variant: {other:?}"),
1028 }
1029
1030 let err = ts
1031 .step(50.0, -1.0)
1032 .expect_err("negative erosion_rate must fail");
1033 match err {
1034 FusionError::ConfigError(msg) => assert!(msg.contains("erosion_rate")),
1035 other => panic!("Unexpected error variant: {other:?}"),
1036 }
1037
1038 ts.dt = 0.0;
1039 let err = ts
1040 .evolve_profiles(20.0)
1041 .expect_err("non-positive dt must fail");
1042 match err {
1043 FusionError::ConfigError(msg) => assert!(msg.contains("dt")),
1044 other => panic!("Unexpected error variant: {other:?}"),
1045 }
1046 }
1047
1048 fn iter_neoclassical_params(n: usize) -> NeoclassicalParams {
1053 let rho = Array1::linspace(0.0, 1.0, n);
1054 let q_profile = Array1::from_shape_fn(n, |i| {
1055 let r = rho[i];
1056 1.0 + 2.0 * r * r });
1058 NeoclassicalParams {
1059 r_major: 6.2,
1060 a_minor: 2.0,
1061 b_toroidal: 5.3,
1062 a_ion: 2.0,
1063 z_eff: 1.5,
1064 q_profile,
1065 }
1066 }
1067
1068 #[test]
1069 fn test_neoclassical_chi_positive() {
1070 let nc = iter_neoclassical_params(50);
1071 for i in 1..50 {
1072 let rho = i as f64 / 49.0;
1073 let chi = chang_hinton_chi(rho, 10.0, 10.0, 1.0 + 2.0 * rho * rho, &nc);
1074 assert!(chi > 0.0, "chi must be > 0 at rho={rho}, got {chi}");
1075 assert!(chi.is_finite(), "chi must be finite at rho={rho}");
1076 }
1077 }
1078
1079 #[test]
1080 fn test_neoclassical_chi_increases_with_q() {
1081 let nc = iter_neoclassical_params(50);
1082 let chi_low_q = chang_hinton_chi(0.5, 10.0, 10.0, 1.0, &nc);
1083 let chi_high_q = chang_hinton_chi(0.5, 10.0, 10.0, 3.0, &nc);
1084 assert!(
1085 chi_high_q > chi_low_q,
1086 "Higher q should give larger chi: q=1→{chi_low_q}, q=3→{chi_high_q}"
1087 );
1088 }
1089
1090 #[test]
1091 fn test_neoclassical_chi_floor_for_invalid() {
1092 let nc = iter_neoclassical_params(50);
1093 let chi = chang_hinton_chi(0.0, 10.0, 10.0, 1.0, &nc);
1095 assert!((chi - CHI_NC_FLOOR).abs() < 1e-12);
1096 let chi = chang_hinton_chi(0.5, -1.0, 10.0, 1.5, &nc);
1098 assert!((chi - CHI_NC_FLOOR).abs() < 1e-12);
1099 }
1100
1101 #[test]
1102 fn test_neoclassical_step_no_panic() {
1103 let mut ts = TransportSolver::new();
1104 let nc = iter_neoclassical_params(50);
1105 ts.set_neoclassical(nc).expect("valid neoclassical config");
1106 ts.step(50.0, 1e14)
1107 .expect("valid transport step with neoclassical");
1108 assert!(ts.profiles.te.iter().all(|v| v.is_finite()));
1109 }
1110
1111 #[test]
1112 fn test_neoclassical_differs_from_constant() {
1113 let mut ts_const = TransportSolver::new();
1114 ts_const.update_transport_model(20.0).unwrap();
1115 let chi_const = ts_const.chi.clone();
1116
1117 let mut ts_nc = TransportSolver::new();
1118 let nc = iter_neoclassical_params(50);
1119 ts_nc.set_neoclassical(nc).unwrap();
1120 ts_nc.update_transport_model(20.0).unwrap();
1121
1122 let max_diff = chi_const
1123 .iter()
1124 .zip(ts_nc.chi.iter())
1125 .map(|(a, b)| (a - b).abs())
1126 .fold(0.0f64, f64::max);
1127 assert!(
1128 max_diff > 0.0,
1129 "Neoclassical chi should differ from constant base"
1130 );
1131 }
1132
1133 #[test]
1134 fn test_neoclassical_rejects_invalid_params() {
1135 let mut ts = TransportSolver::new();
1136 let mut nc = iter_neoclassical_params(50);
1137 nc.r_major = -1.0;
1138 assert!(ts.set_neoclassical(nc).is_err());
1139
1140 let mut nc2 = iter_neoclassical_params(50);
1141 nc2.q_profile = Array1::zeros(10); assert!(ts.set_neoclassical(nc2).is_err());
1143 }
1144
1145 #[test]
1146 fn test_chang_hinton_profile_vectorized_matches_scalar() {
1147 let n = 50;
1148 let rho = Array1::linspace(0.02, 1.0, n);
1149 let t_i = Array1::from_shape_fn(n, |i| {
1150 let r = rho[i];
1151 12.0_f64 * (1.0_f64 - r * r).max(0.05_f64)
1152 });
1153 let n_e = Array1::from_shape_fn(n, |i| 8.0 * (1.0 - 0.6 * rho[i] * rho[i]) + 0.5);
1154 let q = Array1::from_shape_fn(n, |i| 1.0 + 2.5 * rho[i] * rho[i]);
1155 let params = NeoclassicalParams {
1156 q_profile: q.clone(),
1157 ..NeoclassicalParams::default()
1158 };
1159
1160 let vectorized = chang_hinton_chi_profile(&rho, &t_i, &n_e, &q, ¶ms);
1161 for i in 0..n {
1162 let scalar = chang_hinton_chi(rho[i], t_i[i], n_e[i], q[i], ¶ms);
1163 assert!(
1164 (vectorized[i] - scalar).abs() < 1e-12,
1165 "vectorized and scalar mismatch at i={i}: vec={}, scalar={scalar}",
1166 vectorized[i]
1167 );
1168 }
1169 }
1170
1171 #[test]
1172 fn test_gyro_bohm_profile_positive_and_scaling() {
1173 let rho = Array1::linspace(0.02, 1.0, 32);
1174 let t_base = Array1::from_elem(32, 4.0);
1175 let t_hot = Array1::from_elem(32, 8.0);
1176
1177 let chi_base = gyro_bohm_chi_profile(&rho, &t_base, 5.0, 2.0, 0.1);
1178 let chi_hot = gyro_bohm_chi_profile(&rho, &t_hot, 5.0, 2.0, 0.1);
1179 let chi_high_b = gyro_bohm_chi_profile(&rho, &t_hot, 10.0, 2.0, 0.1);
1180
1181 assert!(chi_base.iter().all(|v| v.is_finite() && *v > 0.0));
1182 let ratio_t = chi_hot[10] / chi_base[10].max(1e-12);
1183 let ratio_b = chi_hot[10] / chi_high_b[10].max(1e-12);
1184
1185 assert!(
1186 ratio_t > 2.6,
1187 "Expected T^(3/2) scaling (>2.6 for 4->8 keV), got {ratio_t}"
1188 );
1189 assert!(
1190 ratio_b > 3.5,
1191 "Expected ~B^-2 scaling (>3.5 for 5T->10T), got {ratio_b}"
1192 );
1193 }
1194
1195 #[test]
1196 fn test_sauter_bootstrap_profile_is_mostly_positive() {
1197 let n = 60;
1198 let rho = Array1::linspace(0.01, 1.0, n);
1199 let t_i = Array1::from_shape_fn(n, |i| {
1200 let r = rho[i];
1201 12.0_f64 * (1.0_f64 - r * r).max(0.05_f64) + 0.5_f64
1202 });
1203 let t_e = Array1::from_shape_fn(n, |i| {
1204 let r = rho[i];
1205 10.0_f64 * (1.0_f64 - r * r).max(0.05_f64) + 0.3_f64
1206 });
1207 let n_e = Array1::from_shape_fn(n, |i| {
1208 let r = rho[i];
1209 9.0_f64 * (1.0_f64 - 0.7_f64 * r * r) + 1.0_f64
1210 });
1211 let q = Array1::from_shape_fn(n, |i| {
1212 let r = rho[i];
1213 1.0_f64 + 2.5_f64 * r * r
1214 });
1215 let epsilon = rho.mapv(|r| r * 0.32);
1216
1217 let j_bs = sauter_bootstrap_current_profile(&rho, &t_e, &t_i, &n_e, &q, &epsilon, 5.3);
1218 let positive_fraction =
1219 j_bs.iter().filter(|v| **v > 0.0).count() as f64 / j_bs.len() as f64;
1220 assert!(
1221 positive_fraction > 0.75,
1222 "Expected mostly positive bootstrap current, fraction={positive_fraction:.3}"
1223 );
1224 }
1225
1226 #[test]
1227 fn test_cn_tridiag_diagonal_dominance() {
1228 let n = 50;
1229 let rho = Array1::linspace(0.0, 1.0, n);
1230 let chi = Array1::from_elem(n, 0.5);
1231 let (a, b, c, _) = build_cn_tridiag(&chi, &rho, 0.01);
1232
1233 for i in 1..(n - 1) {
1234 let dominance = b[i].abs() - (a[i].abs() + c[i].abs());
1235 assert!(
1236 dominance > 0.0,
1237 "Expected strict diagonal dominance at i={i}: b={}, a={}, c={}",
1238 b[i],
1239 a[i],
1240 c[i]
1241 );
1242 }
1243 }
1244
1245 fn thermal_energy_proxy(ts: &TransportSolver) -> f64 {
1246 ts.profiles
1247 .rho
1248 .iter()
1249 .zip(ts.profiles.ne.iter())
1250 .zip(ts.profiles.ti.iter())
1251 .map(|((r, n_e), t_i)| r.max(0.0) * n_e.max(0.0) * t_i.max(0.0))
1252 .sum::<f64>()
1253 }
1254
1255 #[test]
1256 fn test_transport_step_conserves_energy_proxy() {
1257 let mut ts = TransportSolver::new();
1258 let p_aux = 20.0;
1259 let dt = 0.01;
1260 let before = thermal_energy_proxy(&ts);
1261 transport_step(&mut ts, p_aux, dt).expect("valid transport step");
1262 let after = thermal_energy_proxy(&ts);
1263
1264 let delta = (after - before).abs();
1265 let bound = p_aux * dt * 5e4;
1267 assert!(
1268 delta <= bound,
1269 "Energy proxy change too large: delta={delta}, bound={bound}"
1270 );
1271 }
1272}