Skip to main content

oxiphysics_materials/
simd_paths.rs

1// Copyright 2026 COOLJAPAN OU (Team KitaSan)
2// SPDX-License-Identifier: Apache-2.0
3
4//! SIMD-accelerated bulk material evaluation using Structure-of-Arrays (SoA) layout.
5//!
6//! This module provides batch (N-point) evaluation kernels for the most
7//! compute-intensive constitutive relations in OxiPhysics.  Data is arranged in
8//! SoA format so that LLVM's auto-vectorizer can emit wide SIMD instructions
9//! (SSE2 / AVX2 / NEON) on the inner loops.
10//!
11//! ## Constitutive kernels
12//!
13//! | Function | Model | Description |
14//! |---|---|---|
15//! | [`elastic_stress_batch`] | Isotropic linear elastic | σ = λ tr(ε) I + 2μ ε |
16//! | [`von_mises_yield_batch`] | Von Mises plasticity | f = √(3/2 s:s) − σ_y |
17//! | [`neo_hookean_stress_batch`] | Neo-Hookean hyperelastic | First Piola–Kirchhoff via F |
18//! | [`viscoplastic_rate_batch`] | Perzyna viscoplasticity | γ̇ = ⟨f/σ_y⟩ⁿ / η |
19//! | [`miner_damage_batch`] | Miner's rule fatigue | D += n/N_f per cycle block |
20//! | [`thermal_expansion_stress_batch`] | Thermo-elastic | σ_th = −3K α ΔT I |
21//!
22//! ## SoA storage
23//!
24//! [`SoaMaterialPoints`] stores the full state of N integration points with each
25//! field in its own contiguous `Vec<f64>` so the batch kernels read and write
26//! aligned slices that the compiler can vectorize.
27//!
28//! ## Usage
29//!
30//! ```
31//! use oxiphysics_materials::simd_paths::{SoaMaterialPoints, elastic_stress_batch};
32//!
33//! let n = 8;
34//! let mut pts = SoaMaterialPoints::new(n);
35//! // Set strains: ε_xx = 0.001 for all points
36//! for i in 0..n { pts.strain_xx[i] = 0.001; }
37//!
38//! let e  = 210e9_f64;
39//! let nu = 0.3;
40//! elastic_stress_batch(&mut pts, e, nu);
41//!
42//! // σ_xx should equal E/(1-ν²) * ε_xx for plane stress
43//! assert!(pts.stress_xx[0] > 0.0);
44//! ```
45
46#![allow(dead_code)]
47
48// ── SoaMaterialPoints ─────────────────────────────────────────────────────────
49
50/// Structure-of-Arrays storage for `n` material integration points.
51///
52/// All arrays have the same length `n` so batch kernels can iterate over them
53/// with simple index loops that LLVM can auto-vectorize.
54#[derive(Debug, Clone, Default)]
55pub struct SoaMaterialPoints {
56    /// Number of material points.
57    pub n: usize,
58
59    // ── Strain tensor components (Voigt notation, symmetric) ──────────────
60    /// Normal strain εxx.
61    pub strain_xx: Vec<f64>,
62    /// Normal strain εyy.
63    pub strain_yy: Vec<f64>,
64    /// Normal strain εzz.
65    pub strain_zz: Vec<f64>,
66    /// Shear strain εxy (engineering = 2 * tensorial).
67    pub strain_xy: Vec<f64>,
68    /// Shear strain εyz.
69    pub strain_yz: Vec<f64>,
70    /// Shear strain εxz.
71    pub strain_xz: Vec<f64>,
72
73    // ── Stress tensor components ──────────────────────────────────────────
74    /// Normal stress σxx (output from constitutive kernels).
75    pub stress_xx: Vec<f64>,
76    /// Normal stress σyy.
77    pub stress_yy: Vec<f64>,
78    /// Normal stress σzz.
79    pub stress_zz: Vec<f64>,
80    /// Shear stress σxy.
81    pub stress_xy: Vec<f64>,
82    /// Shear stress σyz.
83    pub stress_yz: Vec<f64>,
84    /// Shear stress σxz.
85    pub stress_xz: Vec<f64>,
86
87    // ── Scalar state fields ───────────────────────────────────────────────
88    /// Accumulated equivalent plastic strain p̄.
89    pub plastic_strain: Vec<f64>,
90    /// Current yield stress σ_y (may harden).
91    pub yield_stress: Vec<f64>,
92    /// Scalar damage variable D ∈ [0, 1].
93    pub damage: Vec<f64>,
94    /// Temperature (K).
95    pub temperature: Vec<f64>,
96    /// Accumulated fatigue damage (Miner's rule) ∈ [0, 1].
97    pub fatigue_damage: Vec<f64>,
98
99    // ── Deformation gradient (row-major 3×3) ─────────────────────────────
100    /// F₁₁ components (one per point).
101    pub f11: Vec<f64>,
102    /// F₁₂.
103    pub f12: Vec<f64>,
104    /// F₁₃.
105    pub f13: Vec<f64>,
106    /// F₂₁.
107    pub f21: Vec<f64>,
108    /// F₂₂.
109    pub f22: Vec<f64>,
110    /// F₂₃.
111    pub f23: Vec<f64>,
112    /// F₃₁.
113    pub f31: Vec<f64>,
114    /// F₃₂.
115    pub f32: Vec<f64>,
116    /// F₃₃.
117    pub f33: Vec<f64>,
118}
119
120impl SoaMaterialPoints {
121    /// Allocate `n` material points, all zero-initialised.
122    pub fn new(n: usize) -> Self {
123        // Identity deformation gradients
124        let mut pts = Self {
125            n,
126            strain_xx: vec![0.0; n],
127            strain_yy: vec![0.0; n],
128            strain_zz: vec![0.0; n],
129            strain_xy: vec![0.0; n],
130            strain_yz: vec![0.0; n],
131            strain_xz: vec![0.0; n],
132            stress_xx: vec![0.0; n],
133            stress_yy: vec![0.0; n],
134            stress_zz: vec![0.0; n],
135            stress_xy: vec![0.0; n],
136            stress_yz: vec![0.0; n],
137            stress_xz: vec![0.0; n],
138            plastic_strain: vec![0.0; n],
139            yield_stress: vec![250e6; n], // default 250 MPa
140            damage: vec![0.0; n],
141            temperature: vec![293.15; n], // room temperature K
142            fatigue_damage: vec![0.0; n],
143            f11: vec![1.0; n],
144            f12: vec![0.0; n],
145            f13: vec![0.0; n],
146            f21: vec![0.0; n],
147            f22: vec![1.0; n],
148            f23: vec![0.0; n],
149            f31: vec![0.0; n],
150            f32: vec![0.0; n],
151            f33: vec![1.0; n],
152        };
153        // yield_stress default overrides
154        pts.yield_stress.fill(250e6);
155        pts
156    }
157
158    /// Reset all stresses to zero.
159    pub fn zero_stress(&mut self) {
160        self.stress_xx.fill(0.0);
161        self.stress_yy.fill(0.0);
162        self.stress_zz.fill(0.0);
163        self.stress_xy.fill(0.0);
164        self.stress_yz.fill(0.0);
165        self.stress_xz.fill(0.0);
166    }
167
168    /// Compute Von Mises equivalent stress for each point.
169    ///
170    /// σ_vm = √( ½ [(σxx−σyy)² + (σyy−σzz)² + (σzz−σxx)²] + 3(σxy² + σyz² + σxz²) )
171    pub fn von_mises_stress(&self) -> Vec<f64> {
172        let n = self.n;
173        let mut vm = Vec::with_capacity(n);
174        for i in 0..n {
175            let dxy = self.stress_xx[i] - self.stress_yy[i];
176            let dyz = self.stress_yy[i] - self.stress_zz[i];
177            let dzx = self.stress_zz[i] - self.stress_xx[i];
178            let shear = self.stress_xy[i] * self.stress_xy[i]
179                + self.stress_yz[i] * self.stress_yz[i]
180                + self.stress_xz[i] * self.stress_xz[i];
181            vm.push((0.5 * (dxy * dxy + dyz * dyz + dzx * dzx) + 3.0 * shear).sqrt());
182        }
183        vm
184    }
185
186    /// Mean (hydrostatic) stress for each point.
187    pub fn hydrostatic_stress(&self) -> Vec<f64> {
188        (0..self.n)
189            .map(|i| (self.stress_xx[i] + self.stress_yy[i] + self.stress_zz[i]) / 3.0)
190            .collect()
191    }
192}
193
194// ── Lamé parameters ───────────────────────────────────────────────────────────
195
196/// Compute Lamé parameters (λ, μ) from Young's modulus E and Poisson ratio ν.
197#[inline]
198pub fn lame(e: f64, nu: f64) -> (f64, f64) {
199    let mu = e / (2.0 * (1.0 + nu));
200    let lambda = e * nu / ((1.0 + nu) * (1.0 - 2.0 * nu));
201    (lambda, mu)
202}
203
204// ── elastic_stress_batch ──────────────────────────────────────────────────────
205
206/// Compute isotropic linear elastic stress for all `n` material points.
207///
208/// σ = λ tr(ε) I + 2μ ε
209///
210/// Updates `pts.stress_{xx,yy,zz,xy,yz,xz}` in place.  The loop is written in
211/// a form that LLVM can auto-vectorize with 4-wide f64 SIMD lanes.
212pub fn elastic_stress_batch(pts: &mut SoaMaterialPoints, e: f64, nu: f64) {
213    let (lambda, mu) = lame(e, nu);
214    let n = pts.n;
215    let two_mu = 2.0 * mu;
216
217    // Batch loop: auto-vectorized by LLVM into AVX2 (4×f64) or SSE2 (2×f64)
218    for i in 0..n {
219        let tr = pts.strain_xx[i] + pts.strain_yy[i] + pts.strain_zz[i];
220        let ltr = lambda * tr;
221        pts.stress_xx[i] = ltr + two_mu * pts.strain_xx[i];
222        pts.stress_yy[i] = ltr + two_mu * pts.strain_yy[i];
223        pts.stress_zz[i] = ltr + two_mu * pts.strain_zz[i];
224        pts.stress_xy[i] = mu * pts.strain_xy[i]; // engineering shear
225        pts.stress_yz[i] = mu * pts.strain_yz[i];
226        pts.stress_xz[i] = mu * pts.strain_xz[i];
227    }
228}
229
230/// Same as [`elastic_stress_batch`] but processes 4 points at once with
231/// manually unrolled 4-wide scalar arithmetic (auto-vectorizer helper).
232///
233/// Use for very large point clouds where the batch size is always a multiple of 4.
234#[inline(always)]
235pub fn elastic_stress_batch_x4(pts: &mut SoaMaterialPoints, e: f64, nu: f64) {
236    let (lambda, mu) = lame(e, nu);
237    let two_mu = 2.0 * mu;
238    let n = pts.n;
239    let chunks = n / 4;
240
241    for c in 0..chunks {
242        let base = c * 4;
243        // Process 4 points per iteration — LLVM fuses into a single AVX2 kernel
244        let tr0 = pts.strain_xx[base] + pts.strain_yy[base] + pts.strain_zz[base];
245        let tr1 = pts.strain_xx[base + 1] + pts.strain_yy[base + 1] + pts.strain_zz[base + 1];
246        let tr2 = pts.strain_xx[base + 2] + pts.strain_yy[base + 2] + pts.strain_zz[base + 2];
247        let tr3 = pts.strain_xx[base + 3] + pts.strain_yy[base + 3] + pts.strain_zz[base + 3];
248
249        let ltr0 = lambda * tr0;
250        let ltr1 = lambda * tr1;
251        let ltr2 = lambda * tr2;
252        let ltr3 = lambda * tr3;
253
254        pts.stress_xx[base] = ltr0 + two_mu * pts.strain_xx[base];
255        pts.stress_xx[base + 1] = ltr1 + two_mu * pts.strain_xx[base + 1];
256        pts.stress_xx[base + 2] = ltr2 + two_mu * pts.strain_xx[base + 2];
257        pts.stress_xx[base + 3] = ltr3 + two_mu * pts.strain_xx[base + 3];
258
259        pts.stress_yy[base] = ltr0 + two_mu * pts.strain_yy[base];
260        pts.stress_yy[base + 1] = ltr1 + two_mu * pts.strain_yy[base + 1];
261        pts.stress_yy[base + 2] = ltr2 + two_mu * pts.strain_yy[base + 2];
262        pts.stress_yy[base + 3] = ltr3 + two_mu * pts.strain_yy[base + 3];
263
264        pts.stress_zz[base] = ltr0 + two_mu * pts.strain_zz[base];
265        pts.stress_zz[base + 1] = ltr1 + two_mu * pts.strain_zz[base + 1];
266        pts.stress_zz[base + 2] = ltr2 + two_mu * pts.strain_zz[base + 2];
267        pts.stress_zz[base + 3] = ltr3 + two_mu * pts.strain_zz[base + 3];
268
269        pts.stress_xy[base] = mu * pts.strain_xy[base];
270        pts.stress_xy[base + 1] = mu * pts.strain_xy[base + 1];
271        pts.stress_xy[base + 2] = mu * pts.strain_xy[base + 2];
272        pts.stress_xy[base + 3] = mu * pts.strain_xy[base + 3];
273
274        pts.stress_yz[base] = mu * pts.strain_yz[base];
275        pts.stress_yz[base + 1] = mu * pts.strain_yz[base + 1];
276        pts.stress_yz[base + 2] = mu * pts.strain_yz[base + 2];
277        pts.stress_yz[base + 3] = mu * pts.strain_yz[base + 3];
278
279        pts.stress_xz[base] = mu * pts.strain_xz[base];
280        pts.stress_xz[base + 1] = mu * pts.strain_xz[base + 1];
281        pts.stress_xz[base + 2] = mu * pts.strain_xz[base + 2];
282        pts.stress_xz[base + 3] = mu * pts.strain_xz[base + 3];
283    }
284
285    // Scalar tail
286    for i in (chunks * 4)..n {
287        let tr = pts.strain_xx[i] + pts.strain_yy[i] + pts.strain_zz[i];
288        let ltr = lambda * tr;
289        pts.stress_xx[i] = ltr + two_mu * pts.strain_xx[i];
290        pts.stress_yy[i] = ltr + two_mu * pts.strain_yy[i];
291        pts.stress_zz[i] = ltr + two_mu * pts.strain_zz[i];
292        pts.stress_xy[i] = mu * pts.strain_xy[i];
293        pts.stress_yz[i] = mu * pts.strain_yz[i];
294        pts.stress_xz[i] = mu * pts.strain_xz[i];
295    }
296}
297
298// ── von_mises_yield_batch ─────────────────────────────────────────────────────
299
300/// Evaluate the Von Mises yield function for each material point.
301///
302/// `f_i = σ_vm_i − σ_y_i`
303///
304/// * Positive values indicate plastic flow.
305/// * `pts.stress_*` must be populated before calling (e.g. via [`elastic_stress_batch`]).
306///
307/// Returns a `Vec<f64>` of yield function values.
308pub fn von_mises_yield_batch(pts: &SoaMaterialPoints) -> Vec<f64> {
309    let vm = pts.von_mises_stress();
310    (0..pts.n).map(|i| vm[i] - pts.yield_stress[i]).collect()
311}
312
313// ── neo_hookean_stress_batch ──────────────────────────────────────────────────
314
315/// Compute the Cauchy stress for a Neo-Hookean hyperelastic material at each point.
316///
317/// The Neo-Hookean model:
318/// ```text
319/// Ψ = μ/2 (tr(b) − 3) − μ ln(J) + λ/2 (ln J)²
320/// σ = (1/J) [ μ (b − I) + λ ln(J) I ]
321/// ```
322/// where `b = F Fᵀ` (left Cauchy-Green tensor) and `J = det(F)`.
323///
324/// Updates `pts.stress_*` from the stored deformation gradients `pts.f*`.
325pub fn neo_hookean_stress_batch(pts: &mut SoaMaterialPoints, e: f64, nu: f64) {
326    let (lambda, mu) = lame(e, nu);
327    let n = pts.n;
328
329    for i in 0..n {
330        // Deformation gradient F (row-major)
331        let f = [
332            [pts.f11[i], pts.f12[i], pts.f13[i]],
333            [pts.f21[i], pts.f22[i], pts.f23[i]],
334            [pts.f31[i], pts.f32[i], pts.f33[i]],
335        ];
336
337        // J = det(F)
338        let j = det3(f);
339        if j.abs() < 1e-15 {
340            continue; // degenerate element, skip
341        }
342        let ln_j = j.ln();
343        let j_inv = 1.0 / j;
344
345        // b = F Fᵀ (left Cauchy-Green)
346        let b = mat_mat_t(f, f);
347
348        // σ = (1/J) [ μ (b − I) + λ ln(J) I ]
349        pts.stress_xx[i] = j_inv * (mu * (b[0][0] - 1.0) + lambda * ln_j);
350        pts.stress_yy[i] = j_inv * (mu * (b[1][1] - 1.0) + lambda * ln_j);
351        pts.stress_zz[i] = j_inv * (mu * (b[2][2] - 1.0) + lambda * ln_j);
352        pts.stress_xy[i] = j_inv * mu * b[0][1];
353        pts.stress_yz[i] = j_inv * mu * b[1][2];
354        pts.stress_xz[i] = j_inv * mu * b[0][2];
355    }
356}
357
358// ── viscoplastic_rate_batch ───────────────────────────────────────────────────
359
360/// Compute the viscoplastic strain rate for each point (Perzyna model).
361///
362/// `γ̇_i = ⟨f_i / σ_y_i⟩ⁿ / η`
363///
364/// where `⟨x⟩ = max(0, x)` (Macaulay bracket).
365///
366/// * `yield_values` — yield function values (from [`von_mises_yield_batch`])
367/// * `n_exp`        — rate sensitivity exponent (typical 1–10)
368/// * `eta`          — viscosity parameter (typical 1e-3 – 1e-1)
369///
370/// Returns `Vec<f64>` of plastic strain rates.
371pub fn viscoplastic_rate_batch(
372    pts: &SoaMaterialPoints,
373    yield_values: &[f64],
374    n_exp: f64,
375    eta: f64,
376) -> Vec<f64> {
377    assert_eq!(yield_values.len(), pts.n);
378    let eta_inv = 1.0 / eta.max(1e-300);
379    (0..pts.n)
380        .map(|i| {
381            let f = yield_values[i].max(0.0);
382            let overstress = f / pts.yield_stress[i].max(1e-15);
383            overstress.powf(n_exp) * eta_inv
384        })
385        .collect()
386}
387
388// ── miner_damage_batch ────────────────────────────────────────────────────────
389
390/// Accumulate Miner's rule fatigue damage for a block of `n_applied` loading cycles.
391///
392/// `D_i += n_applied / N_f_i`
393///
394/// where `N_f_i` is the number of cycles to failure for the current stress amplitude
395/// computed via the Basquin relation:
396///
397/// `N_f = (σ_f' / σ_a)^(1/b)`
398///
399/// * `sigma_f_prime` — fatigue strength coefficient
400/// * `b`             — Basquin exponent (typical −0.05 to −0.12)
401/// * `n_applied`     — number of cycles applied this block
402///
403/// Clips `D_i` to `[0, 1]`.  Points where `D_i == 1.0` have failed.
404pub fn miner_damage_batch(pts: &mut SoaMaterialPoints, sigma_f_prime: f64, b: f64, n_applied: f64) {
405    let vm = pts.von_mises_stress();
406    for (i, &sigma_v) in vm.iter().enumerate() {
407        let sigma_a = sigma_v.max(1e-15);
408        // N_f from Basquin: σ_a = σ_f' * (2 N_f)^b  →  N_f = ½ (σ_a / σ_f')^(1/b)
409        // b is negative (typical -0.05 to -0.12), so 1/b is also negative, which means
410        // higher stress amplitude → smaller N_f → more damage per cycle (physically correct).
411        let n_f = 0.5 * (sigma_a / sigma_f_prime.max(1e-15)).powf(1.0 / b);
412        if n_f > 1e-15 {
413            pts.fatigue_damage[i] = (pts.fatigue_damage[i] + n_applied / n_f).min(1.0);
414        }
415    }
416}
417
418// ── thermal_expansion_stress_batch ────────────────────────────────────────────
419
420/// Add the thermal expansion stress increment to existing stress.
421///
422/// `Δσ_th = −3 K α ΔT I`
423///
424/// where `K = E / (3(1 − 2ν))` is the bulk modulus, `α` is the coefficient of
425/// thermal expansion, and `ΔT = T − T_ref`.
426///
427/// * `alpha` — coefficient of thermal expansion (1/K), e.g. 12e-6 for steel
428/// * `t_ref` — stress-free reference temperature (K)
429pub fn thermal_expansion_stress_batch(
430    pts: &mut SoaMaterialPoints,
431    e: f64,
432    nu: f64,
433    alpha: f64,
434    t_ref: f64,
435) {
436    let k_bulk = e / (3.0 * (1.0 - 2.0 * nu));
437    let n = pts.n;
438    for i in 0..n {
439        let dt = pts.temperature[i] - t_ref;
440        let th_stress = -3.0 * k_bulk * alpha * dt;
441        pts.stress_xx[i] += th_stress;
442        pts.stress_yy[i] += th_stress;
443        pts.stress_zz[i] += th_stress;
444    }
445}
446
447// ── return_mapping_batch ──────────────────────────────────────────────────────
448
449/// Radial return mapping for isotropic J2 plasticity (batch, linear isotropic hardening).
450///
451/// Updates both the stress and accumulated plastic strain in place.
452///
453/// Algorithm (per point):
454/// 1. Elastic trial stress (already in `pts.stress_*`)
455/// 2. Compute trial Von Mises: `σ_vm_trial`
456/// 3. If trial < yield: elastic step, do nothing
457/// 4. Else: plastic correction — scale deviatoric back to yield surface
458///
459/// * `e`, `nu`  — elastic constants
460/// * `h`        — isotropic hardening modulus (H = dσ_y / dp̄)
461pub fn return_mapping_batch(pts: &mut SoaMaterialPoints, e: f64, nu: f64, h: f64) {
462    let mu = e / (2.0 * (1.0 + nu));
463    let n = pts.n;
464
465    for i in 0..n {
466        let vm_trial = {
467            let dxy = pts.stress_xx[i] - pts.stress_yy[i];
468            let dyz = pts.stress_yy[i] - pts.stress_zz[i];
469            let dzx = pts.stress_zz[i] - pts.stress_xx[i];
470            (0.5 * (dxy * dxy + dyz * dyz + dzx * dzx)
471                + 3.0
472                    * (pts.stress_xy[i] * pts.stress_xy[i]
473                        + pts.stress_yz[i] * pts.stress_yz[i]
474                        + pts.stress_xz[i] * pts.stress_xz[i]))
475                .sqrt()
476        };
477
478        let sigma_y = pts.yield_stress[i];
479        if vm_trial <= sigma_y + 1e-12 {
480            continue; // elastic
481        }
482
483        // Plastic multiplier Δγ (linear hardening)
484        let dg = (vm_trial - sigma_y) / (3.0 * mu + h);
485
486        // Hydrostatic stress (unchanged in J2)
487        let p = (pts.stress_xx[i] + pts.stress_yy[i] + pts.stress_zz[i]) / 3.0;
488
489        // Deviatoric stress components
490        let sx = pts.stress_xx[i] - p;
491        let sy = pts.stress_yy[i] - p;
492        let sz = pts.stress_zz[i] - p;
493        let txy = pts.stress_xy[i];
494        let tyz = pts.stress_yz[i];
495        let txz = pts.stress_xz[i];
496
497        // Scale deviatoric back: s_corrected = s_trial * (1 - 3μΔγ/σ_vm_trial)
498        let scale = (1.0 - 3.0 * mu * dg / vm_trial.max(1e-15)).max(0.0);
499        pts.stress_xx[i] = p + sx * scale;
500        pts.stress_yy[i] = p + sy * scale;
501        pts.stress_zz[i] = p + sz * scale;
502        pts.stress_xy[i] = txy * scale;
503        pts.stress_yz[i] = tyz * scale;
504        pts.stress_xz[i] = txz * scale;
505
506        // Update internal variables
507        pts.plastic_strain[i] += dg;
508        pts.yield_stress[i] = sigma_y + h * dg; // linear isotropic hardening
509    }
510}
511
512// ── Matrix helpers ─────────────────────────────────────────────────────────────
513
514/// Determinant of a 3×3 matrix (row-major).
515#[inline(always)]
516fn det3(m: [[f64; 3]; 3]) -> f64 {
517    m[0][0] * (m[1][1] * m[2][2] - m[1][2] * m[2][1])
518        - m[0][1] * (m[1][0] * m[2][2] - m[1][2] * m[2][0])
519        + m[0][2] * (m[1][0] * m[2][1] - m[1][1] * m[2][0])
520}
521
522/// Matrix product `A * Bᵀ` for 3×3 matrices.
523#[inline(always)]
524fn mat_mat_t(a: [[f64; 3]; 3], b: [[f64; 3]; 3]) -> [[f64; 3]; 3] {
525    let mut c = [[0.0f64; 3]; 3];
526    for i in 0..3 {
527        for j in 0..3 {
528            for k in 0..3 {
529                c[i][j] += a[i][k] * b[j][k]; // A * B^T
530            }
531        }
532    }
533    c
534}
535
536// ── MaterialBatchStats ────────────────────────────────────────────────────────
537
538/// Summary statistics for a batch of material points.
539#[derive(Debug, Clone)]
540pub struct MaterialBatchStats {
541    /// Number of points.
542    pub n: usize,
543    /// Number of plastically active points (yield function > 0).
544    pub plastic_count: usize,
545    /// Number of failed points (damage >= 1.0 or fatigue >= 1.0).
546    pub failed_count: usize,
547    /// Maximum Von Mises stress.
548    pub max_von_mises: f64,
549    /// Mean Von Mises stress.
550    pub mean_von_mises: f64,
551    /// Maximum accumulated plastic strain.
552    pub max_plastic_strain: f64,
553    /// Maximum fatigue damage.
554    pub max_fatigue_damage: f64,
555}
556
557impl MaterialBatchStats {
558    /// Compute statistics from the current state of `pts`.
559    pub fn compute(pts: &SoaMaterialPoints, yield_values: &[f64]) -> Self {
560        let vm = pts.von_mises_stress();
561        let n = pts.n;
562        let max_vm = vm.iter().cloned().fold(0.0_f64, f64::max);
563        let mean_vm = vm.iter().sum::<f64>() / n.max(1) as f64;
564        let plastic_count = yield_values.iter().filter(|&&f| f > 0.0).count();
565        let failed_count = (0..n)
566            .filter(|&i| pts.damage[i] >= 1.0 || pts.fatigue_damage[i] >= 1.0)
567            .count();
568        let max_ps = pts.plastic_strain.iter().cloned().fold(0.0_f64, f64::max);
569        let max_fd = pts.fatigue_damage.iter().cloned().fold(0.0_f64, f64::max);
570        Self {
571            n,
572            plastic_count,
573            failed_count,
574            max_von_mises: max_vm,
575            mean_von_mises: mean_vm,
576            max_plastic_strain: max_ps,
577            max_fatigue_damage: max_fd,
578        }
579    }
580}
581
582// ── Tests ─────────────────────────────────────────────────────────────────────
583
584#[cfg(test)]
585mod tests {
586    use super::*;
587
588    const E_STEEL: f64 = 210e9;
589    const NU_STEEL: f64 = 0.3;
590
591    #[test]
592    fn elastic_uniaxial_stress() {
593        let n = 4;
594        let mut pts = SoaMaterialPoints::new(n);
595        let eps = 0.001;
596        for i in 0..n {
597            pts.strain_xx[i] = eps;
598        }
599        elastic_stress_batch(&mut pts, E_STEEL, NU_STEEL);
600        // σxx = E ε / (1 - ν²) is NOT the 3D formula; in 3D: σxx = λε + 2με = (λ + 2μ)ε
601        let (lambda, mu) = lame(E_STEEL, NU_STEEL);
602        let expected_xx = (lambda + 2.0 * mu) * eps;
603        for i in 0..n {
604            assert!(
605                (pts.stress_xx[i] - expected_xx).abs() < 1.0, // within 1 Pa
606                "σxx[{i}] = {:.3e}, expected {:.3e}",
607                pts.stress_xx[i],
608                expected_xx
609            );
610        }
611    }
612
613    #[test]
614    fn x4_matches_scalar() {
615        let n = 8;
616        let mut pts1 = SoaMaterialPoints::new(n);
617        let mut pts2 = SoaMaterialPoints::new(n);
618        for i in 0..n {
619            pts1.strain_xx[i] = 0.001 * (i + 1) as f64;
620            pts2.strain_xx[i] = pts1.strain_xx[i];
621        }
622        elastic_stress_batch(&mut pts1, E_STEEL, NU_STEEL);
623        elastic_stress_batch_x4(&mut pts2, E_STEEL, NU_STEEL);
624        for i in 0..n {
625            assert!(
626                (pts1.stress_xx[i] - pts2.stress_xx[i]).abs() < 1e-6,
627                "x4 mismatch at {i}: {} vs {}",
628                pts1.stress_xx[i],
629                pts2.stress_xx[i]
630            );
631        }
632    }
633
634    #[test]
635    fn neo_hookean_identity_gives_zero_stress() {
636        let n = 2;
637        let mut pts = SoaMaterialPoints::new(n); // F = I by default
638        neo_hookean_stress_batch(&mut pts, E_STEEL, NU_STEEL);
639        for i in 0..n {
640            assert!(
641                pts.stress_xx[i].abs() < 1e-3,
642                "σxx should be ~0 at identity F"
643            );
644        }
645    }
646
647    #[test]
648    fn return_mapping_elastic_stays_elastic() {
649        let mut pts = SoaMaterialPoints::new(1);
650        pts.strain_xx[0] = 1e-6; // tiny strain, well below yield
651        elastic_stress_batch(&mut pts, E_STEEL, NU_STEEL);
652        let before = pts.plastic_strain[0];
653        return_mapping_batch(&mut pts, E_STEEL, NU_STEEL, 2e9);
654        assert_eq!(pts.plastic_strain[0], before, "should remain elastic");
655    }
656
657    #[test]
658    fn return_mapping_plastic_updates_state() {
659        let mut pts = SoaMaterialPoints::new(1);
660        // Force high stress: set stress directly above yield
661        pts.stress_xx[0] = 500e6; // 500 MPa >> 250 MPa yield
662        pts.yield_stress[0] = 250e6;
663        return_mapping_batch(&mut pts, E_STEEL, NU_STEEL, 0.0); // H=0 (perfect plasticity)
664        // After return mapping, Von Mises should be ≤ yield + small tolerance
665        let vm = pts.von_mises_stress()[0];
666        assert!(
667            vm <= pts.yield_stress[0] + 1e6,
668            "stress not returned to yield surface: σ_vm={vm:.3e}"
669        );
670        assert!(
671            pts.plastic_strain[0] > 0.0,
672            "plastic strain should accumulate"
673        );
674    }
675
676    #[test]
677    fn miner_damage_accumulates() {
678        let mut pts = SoaMaterialPoints::new(2);
679        pts.stress_xx[0] = 300e6; // high amplitude
680        pts.stress_xx[1] = 100e6; // low amplitude
681        miner_damage_batch(&mut pts, 1000e6, -0.1, 100.0);
682        assert!(
683            pts.fatigue_damage[0] > pts.fatigue_damage[1],
684            "high-stress point should accumulate more damage"
685        );
686    }
687
688    #[test]
689    fn thermal_expansion_compressive_on_heating() {
690        let mut pts = SoaMaterialPoints::new(1);
691        // Start with zero stress, heat the point
692        pts.temperature[0] = 393.15; // 120 °C
693        thermal_expansion_stress_batch(&mut pts, E_STEEL, NU_STEEL, 12e-6, 293.15);
694        // Constrained heating → compressive hydrostatic stress
695        assert!(
696            pts.stress_xx[0] < 0.0,
697            "constrained thermal expansion should be compressive"
698        );
699    }
700
701    #[test]
702    fn von_mises_zero_for_hydrostatic() {
703        let mut pts = SoaMaterialPoints::new(1);
704        // Pure hydrostatic: σxx = σyy = σzz = 100 MPa
705        pts.stress_xx[0] = 100e6;
706        pts.stress_yy[0] = 100e6;
707        pts.stress_zz[0] = 100e6;
708        let vm = pts.von_mises_stress();
709        assert!(
710            vm[0].abs() < 1e-6,
711            "hydrostatic stress has zero Von Mises: {}",
712            vm[0]
713        );
714    }
715
716    #[test]
717    fn viscoplastic_rate_zero_in_elastic_domain() {
718        let pts = SoaMaterialPoints::new(2);
719        let yield_values = vec![-10.0, -1.0]; // below yield
720        let rates = viscoplastic_rate_batch(&pts, &yield_values, 3.0, 1e-3);
721        for &r in &rates {
722            assert_eq!(r, 0.0, "viscoplastic rate should be zero below yield");
723        }
724    }
725
726    #[test]
727    fn batch_stats_counts_failed_points() {
728        let mut pts = SoaMaterialPoints::new(3);
729        pts.fatigue_damage[0] = 1.0; // failed
730        pts.fatigue_damage[1] = 0.5; // partial
731        let yv = vec![0.0; 3];
732        let stats = MaterialBatchStats::compute(&pts, &yv);
733        assert_eq!(stats.failed_count, 1);
734        assert_eq!(stats.n, 3);
735    }
736}