Skip to main content

lie_groups/
sun.rs

1//! Generic SU(N) - Special unitary N×N matrices
2//!
3//! This module provides a compile-time generic implementation of SU(N) for arbitrary N.
4//! It elegantly generalizes SU(2) and SU(3) while maintaining type safety and efficiency.
5//!
6//! # Mathematical Structure
7//!
8//! ```text
9//! SU(N) = { U ∈ ℂᴺˣᴺ | U† U = I, det(U) = 1 }
10//! ```
11//!
12//! # Lie Algebra
13//!
14//! The Lie algebra su(N) consists of N×N traceless anti-Hermitian matrices:
15//! ```text
16//! su(N) = { X ∈ ℂᴺˣᴺ | X† = -X, Tr(X) = 0 }
17//! dim(su(N)) = N² - 1
18//! ```
19//!
20//! # Design Philosophy
21//!
22//! - **Type Safety**: Const generics ensure dimension errors are caught at compile time
23//! - **Efficiency**: Lazy matrix construction, SIMD-friendly operations
24//! - **Elegance**: Unified interface for all N (including N=2,3)
25//! - **Generality**: Works for arbitrary N ≥ 2
26//!
27//! # Examples
28//!
29//! ```ignore
30//! use lie_groups::sun::SunAlgebra;
31//! use lie_groups::LieAlgebra;
32//!
33//! // SU(4) for grand unified theories
34//! type Su4Algebra = SunAlgebra<4>;
35//! let x = Su4Algebra::zero();
36//! assert_eq!(Su4Algebra::DIM, 15);  // 4² - 1 = 15
37//!
38//! // Type safety: dimensions checked at compile time
39//! let su2 = SunAlgebra::<2>::basis_element(0);  // dim = 3
40//! let su3 = SunAlgebra::<3>::basis_element(0);  // dim = 8
41//! // su2.add(&su3);  // Compile error! Incompatible types
42//! ```
43//!
44//! # Physics Applications
45//!
46//! - **SU(2)**: Weak force, isospin
47//! - **SU(3)**: Strong force (QCD), color charge
48//! - **SU(4)**: Pati-Salam model, flavor symmetry
49//! - **SU(5)**: Georgi-Glashow GUT
50//! - **SU(6)**: Flavor SU(3) × color SU(2)
51//!
52//! # Performance
53//!
54//! - Algebra operations: O(N²) `[optimal]`
55//! - Matrix construction: O(N²) `[lazy, only when needed]`
56//! - Exponential map: O(N³) via scaling-and-squaring
57//! - Memory: (N²-1)·sizeof(f64) bytes for algebra
58
59use crate::traits::{
60    AntiHermitianByConstruction, Compact, LieAlgebra, LieGroup, SemiSimple, Simple,
61    TracelessByConstruction,
62};
63use ndarray::Array2;
64use num_complex::Complex64;
65use std::fmt;
66use std::marker::PhantomData;
67use std::ops::{Add, Mul, MulAssign, Neg, Sub};
68
69/// Lie algebra su(N) - (N²-1)-dimensional space of traceless anti-Hermitian matrices
70///
71/// # Type Parameter
72///
73/// - `N`: Matrix dimension (must be ≥ 2)
74///
75/// # Representation
76///
77/// Elements are stored as (N²-1) real coefficients corresponding to the generalized
78/// Gell-Mann basis. The basis is constructed systematically:
79///
80/// 1. **Symmetric generators** (N(N-1)/2 elements):
81///    - λᵢⱼ with i < j: has 1 at (i,j) and (j,i)
82///
83/// 2. **Antisymmetric generators** (N(N-1)/2 elements):
84///    - λᵢⱼ with i < j: has -i at (i,j) and +i at (j,i)
85///
86/// 3. **Diagonal generators** (N-1 elements):
87///    - λₖ diagonal with first k entries = 1, (k+1)-th entry = -k
88///
89/// This generalizes the Pauli matrices (N=2) and Gell-Mann matrices (N=3).
90///
91/// # Mathematical Properties
92///
93/// - Hermitian generators: λⱼ† = λⱼ
94/// - Traceless: Tr(λⱼ) = 0
95/// - Normalized: Tr(λᵢλⱼ) = 2δᵢⱼ
96/// - Completeness: {λⱼ/√2} form orthonormal basis for traceless Hermitian matrices
97///
98/// # Memory Layout
99///
100/// For SU(N), we store (N²-1) f64 values in a heap-allocated Vec for N > 4,
101/// or stack-allocated array for N ≤ 4 (common cases).
102#[derive(Clone, Debug, PartialEq)]
103pub struct SunAlgebra<const N: usize> {
104    /// Coefficients in generalized Gell-Mann basis
105    /// Length: N² - 1
106    pub(crate) coefficients: Vec<f64>,
107    _phantom: PhantomData<[(); N]>,
108}
109
110impl<const N: usize> SunAlgebra<N> {
111    /// Dimension of su(N) algebra: N² - 1, valid only for N ≥ 2
112    const DIM: usize = {
113        assert!(
114            N >= 2,
115            "SU(N) requires N >= 2: SU(1) is trivial, SU(0) is undefined"
116        );
117        N * N - 1
118    };
119
120    /// Create new algebra element from coefficients
121    ///
122    /// # Panics
123    ///
124    /// Panics if `coefficients.len() != N² - 1`
125    #[must_use]
126    pub fn new(coefficients: Vec<f64>) -> Self {
127        assert_eq!(
128            coefficients.len(),
129            Self::DIM,
130            "SU({}) algebra requires {} coefficients, got {}",
131            N,
132            Self::DIM,
133            coefficients.len()
134        );
135        Self {
136            coefficients,
137            _phantom: PhantomData,
138        }
139    }
140
141    /// Returns the coefficients in the generalized Gell-Mann basis.
142    #[must_use]
143    pub fn coefficients(&self) -> &[f64] {
144        &self.coefficients
145    }
146
147    /// Convert to N×N anti-Hermitian matrix: X = i·∑ⱼ aⱼ·(λⱼ/2)
148    ///
149    /// This is the fundamental representation in ℂᴺˣᴺ.
150    /// Convention: tr(Tₐ†Tᵦ) = ½δₐᵦ where Tₐ = iλₐ/2.
151    ///
152    /// # Performance
153    ///
154    /// - Time: O(N²)
155    /// - Space: O(N²)
156    /// - Lazy: Only computed when called
157    ///
158    /// # Mathematical Formula
159    ///
160    /// Given coefficients [a₁, ..., a_{N²-1}], returns:
161    /// ```text
162    /// X = i·∑ⱼ aⱼ·(λⱼ/2)
163    /// ```
164    /// where λⱼ are the generalized Gell-Mann matrices with tr(λₐλᵦ) = 2δₐᵦ.
165    #[must_use]
166    pub fn to_matrix(&self) -> Array2<Complex64> {
167        let mut matrix = Array2::zeros((N, N));
168        let i = Complex64::new(0.0, 1.0);
169
170        let mut idx = 0;
171
172        // Symmetric generators: (i,j) with i < j
173        for row in 0..N {
174            for col in (row + 1)..N {
175                let coeff = self.coefficients[idx];
176                matrix[[row, col]] += i * coeff;
177                matrix[[col, row]] += i * coeff;
178                idx += 1;
179            }
180        }
181
182        // Antisymmetric generators: (i,j) with i < j
183        // Standard Gell-Mann: Λ^A_{ij} = -iE_{ij} + iE_{ji}
184        // T = iΛ/2 gives +coeff/2 (real) at (row,col), -coeff/2 at (col,row)
185        for row in 0..N {
186            for col in (row + 1)..N {
187                let coeff = self.coefficients[idx];
188                matrix[[row, col]] += Complex64::new(coeff, 0.0); // +coeff (real)
189                matrix[[col, row]] += Complex64::new(-coeff, 0.0); // -coeff (real)
190                idx += 1;
191            }
192        }
193
194        // Diagonal generators
195        // The k-th diagonal generator (k=0..N-2) has:
196        // - First (k+1) diagonal entries = +1
197        // - (k+2)-th diagonal entry = -(k+1)
198        // - Normalized so Tr(λ²) = 2
199        for k in 0..(N - 1) {
200            let coeff = self.coefficients[idx];
201
202            // Normalization: √(2 / (k+1)(k+2))
203            let k_f = k as f64;
204            let normalization = 2.0 / ((k_f + 1.0) * (k_f + 2.0));
205            let scale = normalization.sqrt();
206
207            // First (k+1) entries: +1
208            for j in 0..=k {
209                matrix[[j, j]] += i * coeff * scale;
210            }
211            // (k+2)-th entry: -(k+1)
212            matrix[[k + 1, k + 1]] += i * coeff * scale * (-(k_f + 1.0));
213
214            idx += 1;
215        }
216
217        // Apply /2 for tr(Tₐ†Tᵦ) = ½δₐᵦ convention
218        matrix.mapv_inplace(|z| z * 0.5);
219        matrix
220    }
221
222    /// Construct algebra element from matrix
223    ///
224    /// Given X ∈ su(N), extract coefficients in Gell-Mann basis.
225    ///
226    /// # Performance
227    ///
228    /// O(N²) time via inner products with basis elements.
229    #[must_use]
230    pub fn from_matrix(matrix: &Array2<Complex64>) -> Self {
231        assert_eq!(matrix.nrows(), N);
232        assert_eq!(matrix.ncols(), N);
233
234        let mut coefficients = vec![0.0; Self::DIM];
235        let mut idx = 0;
236
237        // Convention: X = i·∑ aⱼ·(λⱼ/2), so matrix entries are half
238        // what they would be for the raw Gell-Mann basis.
239        // We extract by reading the matrix and multiplying by 2.
240
241        // Extract symmetric components
242        // λ has 1 at (row,col) and (col,row)
243        // X[row,col] = i·a/2, so a = 2·Im(X[row,col])
244        for row in 0..N {
245            for col in (row + 1)..N {
246                let val = matrix[[row, col]];
247                coefficients[idx] = val.im * 2.0;
248                idx += 1;
249            }
250        }
251
252        // Extract antisymmetric components
253        // Standard Gell-Mann: Λ^A has -i at (row,col), +i at (col,row)
254        // T = iΛ/2 gives +a/2 (real) at (row,col)
255        // so a = 2·Re(X[row,col])
256        for row in 0..N {
257            for col in (row + 1)..N {
258                let val = matrix[[row, col]];
259                coefficients[idx] = val.re * 2.0;
260                idx += 1;
261            }
262        }
263
264        // Extract diagonal components using proper inner product
265        //
266        // The k-th diagonal generator H_k has:
267        //   - entries [[j,j]] = scale_k for j = 0..=k
268        //   - entry [[k+1, k+1]] = -(k+1) * scale_k
269        //   - normalized so Tr(H_k²) = 2
270        //
271        // With /2 convention: X_diag = i·a_k·(H_k/2)
272        // So a_k = Im(Tr(X · H_k)) / (Tr(H_k²)/2) = Im(Tr(X · H_k))
273        for k in 0..(N - 1) {
274            let k_f = k as f64;
275            let normalization = 2.0 / ((k_f + 1.0) * (k_f + 2.0));
276            let scale = normalization.sqrt();
277
278            let mut trace_prod = Complex64::new(0.0, 0.0);
279            for j in 0..=k {
280                trace_prod += matrix[[j, j]] * scale;
281            }
282            trace_prod += matrix[[k + 1, k + 1]] * (-(k_f + 1.0) * scale);
283
284            // With /2 convention, the extraction formula becomes Im(Tr(X·H_k))
285            // since X = i·a·H_k/2 gives Tr(X·H_k) = i·a·Tr(H_k²)/2 = i·a
286            coefficients[idx] = trace_prod.im;
287            idx += 1;
288        }
289
290        Self::new(coefficients)
291    }
292}
293
294impl<const N: usize> Add for SunAlgebra<N> {
295    type Output = Self;
296    fn add(self, rhs: Self) -> Self {
297        let coefficients = self
298            .coefficients
299            .iter()
300            .zip(&rhs.coefficients)
301            .map(|(a, b)| a + b)
302            .collect();
303        Self::new(coefficients)
304    }
305}
306
307impl<const N: usize> Add<&SunAlgebra<N>> for SunAlgebra<N> {
308    type Output = SunAlgebra<N>;
309    fn add(self, rhs: &SunAlgebra<N>) -> SunAlgebra<N> {
310        let coefficients = self
311            .coefficients
312            .iter()
313            .zip(&rhs.coefficients)
314            .map(|(a, b)| a + b)
315            .collect();
316        Self::new(coefficients)
317    }
318}
319
320impl<const N: usize> Add<SunAlgebra<N>> for &SunAlgebra<N> {
321    type Output = SunAlgebra<N>;
322    fn add(self, rhs: SunAlgebra<N>) -> SunAlgebra<N> {
323        let coefficients = self
324            .coefficients
325            .iter()
326            .zip(&rhs.coefficients)
327            .map(|(a, b)| a + b)
328            .collect();
329        SunAlgebra::<N>::new(coefficients)
330    }
331}
332
333impl<const N: usize> Add<&SunAlgebra<N>> for &SunAlgebra<N> {
334    type Output = SunAlgebra<N>;
335    fn add(self, rhs: &SunAlgebra<N>) -> SunAlgebra<N> {
336        let coefficients = self
337            .coefficients
338            .iter()
339            .zip(&rhs.coefficients)
340            .map(|(a, b)| a + b)
341            .collect();
342        SunAlgebra::<N>::new(coefficients)
343    }
344}
345
346impl<const N: usize> Sub for SunAlgebra<N> {
347    type Output = Self;
348    fn sub(self, rhs: Self) -> Self {
349        let coefficients = self
350            .coefficients
351            .iter()
352            .zip(&rhs.coefficients)
353            .map(|(a, b)| a - b)
354            .collect();
355        Self::new(coefficients)
356    }
357}
358
359impl<const N: usize> Neg for SunAlgebra<N> {
360    type Output = Self;
361    fn neg(self) -> Self {
362        let coefficients = self.coefficients.iter().map(|x| -x).collect();
363        Self::new(coefficients)
364    }
365}
366
367impl<const N: usize> Mul<f64> for SunAlgebra<N> {
368    type Output = Self;
369    fn mul(self, scalar: f64) -> Self {
370        let coefficients = self.coefficients.iter().map(|x| x * scalar).collect();
371        Self::new(coefficients)
372    }
373}
374
375impl<const N: usize> Mul<SunAlgebra<N>> for f64 {
376    type Output = SunAlgebra<N>;
377    fn mul(self, rhs: SunAlgebra<N>) -> SunAlgebra<N> {
378        rhs * self
379    }
380}
381
382impl<const N: usize> LieAlgebra for SunAlgebra<N> {
383    // Model-theoretic guard: SU(N) requires N ≥ 2.
384    // SU(1) = {I} is trivial (algebra is zero-dimensional, bracket undefined).
385    // SU(0) underflows usize arithmetic.
386    // This const assert promotes the degenerate-model failure from runtime panic
387    // to a compile-time error, consistent with the sealed-trait philosophy.
388    const DIM: usize = {
389        assert!(
390            N >= 2,
391            "SU(N) requires N >= 2: SU(1) is trivial, SU(0) is undefined"
392        );
393        N * N - 1
394    };
395
396    fn zero() -> Self {
397        Self {
398            coefficients: vec![0.0; Self::DIM],
399            _phantom: PhantomData,
400        }
401    }
402
403    fn add(&self, other: &Self) -> Self {
404        let coefficients = self
405            .coefficients
406            .iter()
407            .zip(&other.coefficients)
408            .map(|(a, b)| a + b)
409            .collect();
410        Self::new(coefficients)
411    }
412
413    fn scale(&self, scalar: f64) -> Self {
414        let coefficients = self.coefficients.iter().map(|x| x * scalar).collect();
415        Self::new(coefficients)
416    }
417
418    fn norm(&self) -> f64 {
419        self.coefficients
420            .iter()
421            .map(|x| x.powi(2))
422            .sum::<f64>()
423            .sqrt()
424    }
425
426    fn basis_element(i: usize) -> Self {
427        assert!(
428            i < Self::DIM,
429            "Basis index {} out of range for SU({})",
430            i,
431            N
432        );
433        let mut coefficients = vec![0.0; Self::DIM];
434        coefficients[i] = 1.0;
435        Self::new(coefficients)
436    }
437
438    fn from_components(components: &[f64]) -> Self {
439        assert_eq!(
440            components.len(),
441            Self::DIM,
442            "Expected {} components for SU({}), got {}",
443            Self::DIM,
444            N,
445            components.len()
446        );
447        Self::new(components.to_vec())
448    }
449
450    fn to_components(&self) -> Vec<f64> {
451        self.coefficients.clone()
452    }
453
454    /// Lie bracket: [X, Y] = XY - YX
455    ///
456    /// Computed via matrix commutator for generality.
457    ///
458    /// # Performance
459    ///
460    /// - Time: O(N³) [matrix multiplication]
461    /// - Space: O(N²)
462    ///
463    /// # Note
464    ///
465    /// For N=2,3, specialized implementations with structure constants
466    /// would be faster (O(1) and O(1) respectively). This generic version
467    /// prioritizes correctness and simplicity.
468    ///
469    /// # Mathematical Formula
470    ///
471    /// ```text
472    /// [X, Y] = XY - YX
473    /// ```
474    ///
475    /// This satisfies:
476    /// - Antisymmetry: `[X,Y] = -[Y,X]`
477    /// - Jacobi identity: `[X,[Y,Z]] + [Y,[Z,X]] + [Z,[X,Y]] = 0`
478    /// - Bilinearity
479    fn bracket(&self, other: &Self) -> Self {
480        let x = self.to_matrix();
481        let y = other.to_matrix();
482        let commutator = x.dot(&y) - y.dot(&x);
483        Self::from_matrix(&commutator)
484    }
485
486    #[inline]
487    fn inner(&self, other: &Self) -> f64 {
488        self.coefficients
489            .iter()
490            .zip(other.coefficients.iter())
491            .map(|(a, b)| a * b)
492            .sum()
493    }
494}
495
496/// SU(N) group element - N×N unitary matrix with det = 1
497///
498/// # Type Parameter
499///
500/// - `N`: Matrix dimension
501///
502/// # Representation
503///
504/// Stored as N×N complex matrix satisfying:
505/// - U†U = I (unitarity)
506/// - det(U) = 1 (special)
507///
508/// # Verification
509///
510/// Use `verify_unitarity()` to check constraints numerically.
511#[derive(Debug, Clone, PartialEq)]
512pub struct SUN<const N: usize> {
513    /// N×N complex unitary matrix
514    pub(crate) matrix: Array2<Complex64>,
515}
516
517impl<const N: usize> SUN<N> {
518    /// Access the underlying N×N unitary matrix
519    #[must_use]
520    pub fn matrix(&self) -> &Array2<Complex64> {
521        &self.matrix
522    }
523
524    /// Identity element: Iₙ
525    #[must_use]
526    pub fn identity() -> Self {
527        Self {
528            matrix: Array2::eye(N),
529        }
530    }
531
532    /// Trace of the matrix: Tr(U) = Σᵢ Uᵢᵢ
533    #[must_use]
534    pub fn trace(&self) -> Complex64 {
535        (0..N).map(|i| self.matrix[[i, i]]).sum()
536    }
537
538    /// Verify unitarity: ||U†U - I|| < ε
539    ///
540    /// # Arguments
541    ///
542    /// - `tolerance`: Maximum Frobenius norm deviation
543    ///
544    /// # Returns
545    ///
546    /// `true` if U†U ≈ I within tolerance
547    #[must_use]
548    pub fn verify_unitarity(&self, tolerance: f64) -> bool {
549        let adjoint = self.matrix.t().mapv(|z| z.conj());
550        let product = adjoint.dot(&self.matrix);
551        let identity: Array2<Complex64> = Array2::eye(N);
552
553        let diff = &product - &identity;
554        let norm_sq: f64 = diff.iter().map(num_complex::Complex::norm_sqr).sum();
555
556        norm_sq.sqrt() < tolerance
557    }
558
559    /// Compute determinant
560    ///
561    /// For SU(N), the determinant should be exactly 1 by definition.
562    ///
563    /// # Implementation
564    ///
565    /// - **N=2**: Direct formula `ad - bc`
566    /// - **N=3**: Sarrus' rule / cofactor expansion
567    /// - **N>3**: Returns `1.0` (assumes matrix is on SU(N) manifold)
568    ///
569    /// # Limitations
570    ///
571    /// For N > 3, this function **does not compute the actual determinant**.
572    /// It returns 1.0 under the assumption that matrices constructed via
573    /// `exp()` or `reorthogonalize()` remain on the SU(N) manifold.
574    ///
575    /// To verify unitarity for N > 3, use `verify_special_unitarity()` instead,
576    /// which checks `U†U = I` (a stronger condition than det=1).
577    ///
578    /// For actual determinant computation with N > 3, enable the `ndarray-linalg`
579    /// feature (not currently available) or compute via eigenvalue product.
580    #[must_use]
581    #[allow(clippy::many_single_char_names)] // Standard math notation for matrix elements
582    pub fn determinant(&self) -> Complex64 {
583        // For small N, compute directly using Leibniz formula
584        if N == 2 {
585            let a = self.matrix[[0, 0]];
586            let b = self.matrix[[0, 1]];
587            let c = self.matrix[[1, 0]];
588            let d = self.matrix[[1, 1]];
589            return a * d - b * c;
590        }
591
592        if N == 3 {
593            // 3x3 determinant via Sarrus' rule / cofactor expansion
594            let (a, b, c, d, e, f, g, h, i) = {
595                let m = &self.matrix;
596                (
597                    m[[0, 0]],
598                    m[[0, 1]],
599                    m[[0, 2]],
600                    m[[1, 0]],
601                    m[[1, 1]],
602                    m[[1, 2]],
603                    m[[2, 0]],
604                    m[[2, 1]],
605                    m[[2, 2]],
606                )
607            };
608
609            // det = a(ei - fh) - b(di - fg) + c(dh - eg)
610            return a * (e * i - f * h) - b * (d * i - f * g) + c * (d * h - e * g);
611        }
612
613        // For N > 3: LU decomposition would be ideal, but requires ndarray-linalg
614        // For now, return 1.0 since matrices constructed via exp() preserve det=1
615        // This is valid for elements on the SU(N) manifold
616        Complex64::new(1.0, 0.0)
617    }
618
619    /// Gram-Schmidt reorthogonalization for SU(N) matrices
620    ///
621    /// Projects a potentially corrupted matrix back onto the SU(N) manifold
622    /// using Gram-Schmidt orthogonalization followed by determinant correction.
623    ///
624    /// # Algorithm
625    ///
626    /// 1. Orthogonalize columns using Modified Gram-Schmidt (MGS)
627    /// 2. Normalize to ensure unitarity
628    /// 3. Adjust phase to ensure det(U) = 1
629    ///
630    /// This avoids the log-exp round-trip that would cause infinite recursion
631    /// when called from within `exp()`.
632    ///
633    /// # Numerical Stability
634    ///
635    /// Uses Modified Gram-Schmidt (not Classical GS) for better numerical stability.
636    /// Projections are computed against already-orthonormalized vectors and
637    /// subtracted immediately, providing O(ε) backward error vs O(κε) for CGS.
638    ///
639    /// Reference: Björck, "Numerical Methods for Least Squares Problems" (1996)
640    #[must_use]
641    fn gram_schmidt_project(matrix: Array2<Complex64>) -> Array2<Complex64> {
642        let mut result: Array2<Complex64> = Array2::zeros((N, N));
643
644        // Modified Gram-Schmidt on columns
645        for j in 0..N {
646            let mut col = matrix.column(j).to_owned();
647
648            // Subtract projections onto previous columns
649            for k in 0..j {
650                let prev_col = result.column(k);
651                let proj: Complex64 = prev_col
652                    .iter()
653                    .zip(col.iter())
654                    .map(|(p, c)| p.conj() * c)
655                    .sum();
656                for i in 0..N {
657                    col[i] -= proj * prev_col[i];
658                }
659            }
660
661            // Normalize
662            let norm: f64 = col
663                .iter()
664                .map(num_complex::Complex::norm_sqr)
665                .sum::<f64>()
666                .sqrt();
667
668            // Detect linear dependence: column became zero after orthogonalization
669            debug_assert!(
670                norm > 1e-14,
671                "Gram-Schmidt: column {} is linearly dependent (norm = {:.2e}). \
672                 Input matrix is rank-deficient.",
673                j,
674                norm
675            );
676
677            if norm > 1e-14 {
678                for i in 0..N {
679                    result[[i, j]] = col[i] / norm;
680                }
681            }
682            // Note: if norm ≤ 1e-14, column remains zero → det will be ~0 → identity fallback
683        }
684
685        // For N=2 or N=3, compute determinant and fix phase
686        // For larger N, approximate (Gram-Schmidt usually produces det ≈ 1 already)
687        if N <= 3 {
688            let det = if N == 2 {
689                result[[0, 0]] * result[[1, 1]] - result[[0, 1]] * result[[1, 0]]
690            } else {
691                // N=3
692                result[[0, 0]] * (result[[1, 1]] * result[[2, 2]] - result[[1, 2]] * result[[2, 1]])
693                    - result[[0, 1]]
694                        * (result[[1, 0]] * result[[2, 2]] - result[[1, 2]] * result[[2, 0]])
695                    + result[[0, 2]]
696                        * (result[[1, 0]] * result[[2, 1]] - result[[1, 1]] * result[[2, 0]])
697            };
698
699            // Guard against zero determinant (degenerate matrix)
700            let det_norm = det.norm();
701            if det_norm < 1e-14 {
702                // Matrix is degenerate; return identity as fallback
703                return Array2::eye(N);
704            }
705
706            let det_phase = det / det_norm;
707            let correction = (det_phase.conj()).powf(1.0 / N as f64);
708            result.mapv_inplace(|z| z * correction);
709        }
710
711        result
712    }
713
714    /// Distance from identity: ||U - I||_F
715    ///
716    /// Frobenius norm of difference from identity.
717    #[must_use]
718    pub fn distance_to_identity(&self) -> f64 {
719        let identity: Array2<Complex64> = Array2::eye(N);
720        let diff = &self.matrix - &identity;
721        diff.iter()
722            .map(num_complex::Complex::norm_sqr)
723            .sum::<f64>()
724            .sqrt()
725    }
726}
727
728impl<const N: usize> approx::AbsDiffEq for SunAlgebra<N> {
729    type Epsilon = f64;
730
731    fn default_epsilon() -> Self::Epsilon {
732        1e-10
733    }
734
735    fn abs_diff_eq(&self, other: &Self, epsilon: Self::Epsilon) -> bool {
736        self.coefficients
737            .iter()
738            .zip(other.coefficients.iter())
739            .all(|(a, b)| (a - b).abs() < epsilon)
740    }
741}
742
743impl<const N: usize> approx::RelativeEq for SunAlgebra<N> {
744    fn default_max_relative() -> Self::Epsilon {
745        1e-10
746    }
747
748    fn relative_eq(
749        &self,
750        other: &Self,
751        epsilon: Self::Epsilon,
752        max_relative: Self::Epsilon,
753    ) -> bool {
754        self.coefficients
755            .iter()
756            .zip(other.coefficients.iter())
757            .all(|(a, b)| approx::RelativeEq::relative_eq(a, b, epsilon, max_relative))
758    }
759}
760
761impl<const N: usize> fmt::Display for SunAlgebra<N> {
762    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
763        write!(f, "su({})[", N)?;
764        for (i, c) in self.coefficients.iter().enumerate() {
765            if i > 0 {
766                write!(f, ", ")?;
767            }
768            write!(f, "{:.4}", c)?;
769        }
770        write!(f, "]")
771    }
772}
773
774impl<const N: usize> fmt::Display for SUN<N> {
775    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
776        let dist = self.distance_to_identity();
777        write!(f, "SU({})(d={:.4})", N, dist)
778    }
779}
780
781/// Group multiplication: U₁ · U₂
782impl<const N: usize> Mul<&SUN<N>> for &SUN<N> {
783    type Output = SUN<N>;
784    fn mul(self, rhs: &SUN<N>) -> SUN<N> {
785        SUN {
786            matrix: self.matrix.dot(&rhs.matrix),
787        }
788    }
789}
790
791impl<const N: usize> Mul<&SUN<N>> for SUN<N> {
792    type Output = SUN<N>;
793    fn mul(self, rhs: &SUN<N>) -> SUN<N> {
794        &self * rhs
795    }
796}
797
798impl<const N: usize> MulAssign<&SUN<N>> for SUN<N> {
799    fn mul_assign(&mut self, rhs: &SUN<N>) {
800        self.matrix = self.matrix.dot(&rhs.matrix);
801    }
802}
803
804impl<const N: usize> LieGroup for SUN<N> {
805    type Algebra = SunAlgebra<N>;
806    const MATRIX_DIM: usize = {
807        assert!(
808            N >= 2,
809            "SU(N) requires N >= 2: SU(1) is trivial, SU(0) is undefined"
810        );
811        N
812    };
813
814    fn identity() -> Self {
815        Self::identity()
816    }
817
818    fn compose(&self, other: &Self) -> Self {
819        Self {
820            matrix: self.matrix.dot(&other.matrix),
821        }
822    }
823
824    fn inverse(&self) -> Self {
825        // For unitary matrices: U⁻¹ = U†
826        Self {
827            matrix: self.matrix.t().mapv(|z| z.conj()),
828        }
829    }
830
831    fn conjugate_transpose(&self) -> Self {
832        Self {
833            matrix: self.matrix.t().mapv(|z| z.conj()),
834        }
835    }
836
837    /// Adjoint action: `Ad_g(X)` = gXg⁻¹
838    ///
839    /// # Mathematical Formula
840    ///
841    /// For g ∈ SU(N) and X ∈ su(N):
842    /// ```text
843    /// Ad_g(X) = gXg⁻¹
844    /// ```
845    ///
846    /// # Properties
847    ///
848    /// - Group homomorphism: Ad_{gh} = `Ad_g` ∘ `Ad_h`
849    /// - Preserves bracket: `Ad_g([X,Y])` = `[Ad_g(X), Ad_g(Y)]`
850    /// - Preserves norm (SU(N) is compact)
851    fn adjoint_action(&self, algebra_element: &Self::Algebra) -> Self::Algebra {
852        let x = algebra_element.to_matrix();
853        let g_inv = self.inverse();
854
855        // Compute gXg⁻¹
856        let result = self.matrix.dot(&x).dot(&g_inv.matrix);
857
858        SunAlgebra::from_matrix(&result)
859    }
860
861    fn distance_to_identity(&self) -> f64 {
862        self.distance_to_identity()
863    }
864
865    /// Exponential map: exp: su(N) → SU(N)
866    ///
867    /// # Algorithm: Scaling-and-Squaring
868    ///
869    /// For X ∈ su(N) with ||X|| large:
870    /// 1. Scale: X' = X / 2^k such that ||X'|| ≤ 0.5
871    /// 2. Taylor: exp(X') ≈ ∑_{n=0}^{15} X'^n / n!
872    /// 3. Square: exp(X) = [exp(X')]^{2^k}
873    ///
874    /// # Properties
875    ///
876    /// - Preserves unitarity: exp(X) ∈ SU(N) for X ∈ su(N)
877    /// - Preserves det = 1 (Tr(X) = 0 ⟹ det(exp(X)) = exp(Tr(X)) = 1)
878    /// - Numerically stable for all ||X||
879    ///
880    /// # Performance
881    ///
882    /// - Time: O(N³·log(||X||))
883    /// - Space: O(N²)
884    ///
885    /// # Accuracy
886    ///
887    /// - Relative error: ~10⁻¹⁴ (double precision)
888    /// - Unitarity preserved to ~10⁻¹²
889    fn exp(tangent: &Self::Algebra) -> Self {
890        let x = tangent.to_matrix();
891        let norm = matrix_frobenius_norm(&x);
892
893        // Determine scaling factor: k such that ||X/2^k|| ≤ 0.5
894        let k = if norm > 0.5 {
895            (norm / 0.5).log2().ceil() as u32
896        } else {
897            0
898        };
899
900        // Scale down
901        let scale_factor = 2_f64.powi(-(k as i32));
902        let x_scaled = x.mapv(|z| z * scale_factor);
903
904        // Taylor series: exp(X') ≈ ∑ X'^n / n!
905        let exp_scaled = matrix_exp_taylor(&x_scaled, 15);
906
907        // Square k times: exp(X) = [exp(X/2^k)]^{2^k}
908        //
909        // Reorthogonalize after EVERY squaring to prevent numerical drift.
910        //
911        // Each matrix multiply accumulates O(nε) orthogonality loss.
912        // After k squarings without correction: O(2^k · nε) total error.
913        // For k=10, n=4: error ~4000ε ≈ 9e-13 — approaching catastrophic loss.
914        //
915        // Cost of Gram-Schmidt is O(N³), same as the matrix multiply itself,
916        // so the overhead factor is ≤2× (negligible vs correctness).
917        let mut result = exp_scaled;
918        for _ in 0..k {
919            result = result.dot(&result);
920            result = Self::gram_schmidt_project(result);
921        }
922
923        Self { matrix: result }
924    }
925
926    fn log(&self) -> crate::error::LogResult<Self::Algebra> {
927        use crate::error::LogError;
928
929        // Matrix logarithm for SU(N) using inverse scaling-squaring algorithm.
930        //
931        // Algorithm (Higham, "Functions of Matrices", Ch. 11):
932        // 1. Take square roots until ||U^{1/2^k} - I|| < 0.5
933        // 2. Use Taylor series for log(I + X) with ||X|| < 0.5 (fast convergence)
934        // 3. Scale back: log(U) = 2^k × log(U^{1/2^k})
935
936        let dist = self.distance_to_identity();
937        const MAX_DISTANCE: f64 = 2.0;
938
939        if dist > MAX_DISTANCE {
940            return Err(LogError::NotNearIdentity {
941                distance: dist,
942                threshold: MAX_DISTANCE,
943            });
944        }
945
946        if dist < 1e-14 {
947            return Ok(SunAlgebra::zero());
948        }
949
950        // Phase 1: Inverse scaling via matrix square roots
951        let identity_matrix: Array2<Complex64> = Array2::eye(N);
952        let mut current = self.matrix.clone();
953        let mut num_sqrts: u32 = 0;
954        const MAX_SQRTS: u32 = 32;
955        const TARGET_NORM: f64 = 0.5;
956
957        while num_sqrts < MAX_SQRTS {
958            let x_matrix = &current - &identity_matrix;
959            let x_norm = matrix_frobenius_norm(&x_matrix);
960            if x_norm < TARGET_NORM {
961                break;
962            }
963            current = matrix_sqrt_db(&current);
964            num_sqrts += 1;
965        }
966
967        // Phase 2: Taylor series for log(I + X) with ||X|| < 0.5
968        let x_matrix = &current - &identity_matrix;
969        let log_matrix = matrix_log_taylor(&x_matrix, 30);
970
971        // Phase 3: Scale back: log(U) = 2^k × log(U^{1/2^k})
972        let scale_factor = (1_u64 << num_sqrts) as f64;
973        let log_scaled = log_matrix.mapv(|z| z * scale_factor);
974
975        Ok(SunAlgebra::from_matrix(&log_scaled))
976    }
977}
978
979/// Compute Frobenius norm of complex matrix: ||A||_F = √(Tr(A†A))
980fn matrix_frobenius_norm(matrix: &Array2<Complex64>) -> f64 {
981    matrix
982        .iter()
983        .map(num_complex::Complex::norm_sqr)
984        .sum::<f64>()
985        .sqrt()
986}
987
988/// Compute matrix exponential via Taylor series
989///
990/// exp(X) = I + X + X²/2! + X³/3! + ... + X^n/n!
991///
992/// Converges rapidly for ||X|| ≤ 0.5.
993fn matrix_exp_taylor(matrix: &Array2<Complex64>, terms: usize) -> Array2<Complex64> {
994    let n = matrix.nrows();
995    let mut result = Array2::eye(n); // I
996    let mut term = Array2::eye(n); // Current term: X^k / k!
997
998    for k in 1..=terms {
999        // term = term · X / k
1000        term = term.dot(matrix).mapv(|z| z / (k as f64));
1001        result += &term;
1002    }
1003
1004    result
1005}
1006
1007/// Compute matrix logarithm via Taylor series
1008///
1009/// log(I + X) = X - X²/2 + X³/3 - X⁴/4 + ... + (-1)^{n+1}·X^n/n
1010///
1011/// Converges for spectral radius ρ(X) < 1. For ||X||_F < 0.5, convergence
1012/// is rapid (30 terms gives ~3e-11 truncation error).
1013fn matrix_log_taylor(matrix: &Array2<Complex64>, terms: usize) -> Array2<Complex64> {
1014    let mut result = matrix.clone(); // First term: X
1015    let mut x_power = matrix.clone(); // Current power: X^k
1016
1017    for k in 2..=terms {
1018        // x_power = X^k
1019        x_power = x_power.dot(matrix);
1020
1021        // Coefficient: (-1)^{k+1} / k
1022        let sign = if k % 2 == 0 { -1.0 } else { 1.0 };
1023        let coefficient = sign / (k as f64);
1024
1025        // Add term to result
1026        result = result + x_power.mapv(|z| z * coefficient);
1027    }
1028
1029    result
1030}
1031
1032/// Complex N×N matrix inverse via Gauss-Jordan elimination with partial pivoting.
1033///
1034/// Returns None if the matrix is singular (pivot < 1e-15).
1035fn matrix_inverse(a: &Array2<Complex64>) -> Option<Array2<Complex64>> {
1036    let n = a.nrows();
1037    assert_eq!(n, a.ncols());
1038
1039    // Build augmented matrix [A | I]
1040    let mut aug = Array2::<Complex64>::zeros((n, 2 * n));
1041    for i in 0..n {
1042        for j in 0..n {
1043            aug[[i, j]] = a[[i, j]];
1044        }
1045        aug[[i, n + i]] = Complex64::new(1.0, 0.0);
1046    }
1047
1048    for col in 0..n {
1049        // Partial pivoting: find row with largest magnitude in this column
1050        let mut max_norm = 0.0;
1051        let mut max_row = col;
1052        for row in col..n {
1053            let norm = aug[[row, col]].norm();
1054            if norm > max_norm {
1055                max_norm = norm;
1056                max_row = row;
1057            }
1058        }
1059        if max_norm < 1e-15 {
1060            return None; // Singular
1061        }
1062
1063        // Swap rows
1064        if max_row != col {
1065            for j in 0..2 * n {
1066                let tmp = aug[[col, j]];
1067                aug[[col, j]] = aug[[max_row, j]];
1068                aug[[max_row, j]] = tmp;
1069            }
1070        }
1071
1072        // Scale pivot row
1073        let pivot = aug[[col, col]];
1074        for j in 0..2 * n {
1075            aug[[col, j]] /= pivot;
1076        }
1077
1078        // Eliminate column in all other rows
1079        for row in 0..n {
1080            if row != col {
1081                let factor = aug[[row, col]];
1082                // Read pivot row values first to avoid borrow conflict
1083                let pivot_row: Vec<Complex64> = (0..2 * n).map(|j| aug[[col, j]]).collect();
1084                for j in 0..2 * n {
1085                    aug[[row, j]] -= factor * pivot_row[j];
1086                }
1087            }
1088        }
1089    }
1090
1091    // Extract inverse from right half
1092    let mut result = Array2::zeros((n, n));
1093    for i in 0..n {
1094        for j in 0..n {
1095            result[[i, j]] = aug[[i, n + j]];
1096        }
1097    }
1098    Some(result)
1099}
1100
1101/// Matrix square root via Denman-Beavers iteration.
1102///
1103/// Computes U^{1/2} for a matrix U close to identity.
1104/// Uses the iteration Y_{k+1} = (Y_k + Z_k^{-1})/2, Z_{k+1} = (Z_k + Y_k^{-1})/2
1105/// which converges quadratically to U^{1/2}.
1106fn matrix_sqrt_db(u: &Array2<Complex64>) -> Array2<Complex64> {
1107    let n = u.nrows();
1108    let mut y = u.clone();
1109    let mut z = Array2::<Complex64>::eye(n);
1110
1111    const MAX_ITERS: usize = 20;
1112    const TOL: f64 = 1e-14;
1113
1114    for _ in 0..MAX_ITERS {
1115        let y_inv = matrix_inverse(&y).unwrap_or_else(|| y.t().mapv(|z| z.conj()));
1116        let z_inv = matrix_inverse(&z).unwrap_or_else(|| z.t().mapv(|z| z.conj()));
1117
1118        let y_new = (&y + &z_inv).mapv(|z| z * 0.5);
1119        let z_new = (&z + &y_inv).mapv(|z| z * 0.5);
1120
1121        let diff = matrix_frobenius_norm(&(&y_new - &y));
1122        y = y_new;
1123        z = z_new;
1124
1125        if diff < TOL {
1126            break;
1127        }
1128    }
1129
1130    y
1131}
1132
1133// ============================================================================
1134// Const Generic Specializations
1135// ============================================================================
1136
1137// ============================================================================
1138// Algebra Marker Traits
1139// ============================================================================
1140
1141/// SU(N) is compact for all N ≥ 2.
1142impl<const N: usize> Compact for SUN<N> {}
1143
1144/// SU(N) is simple for all N ≥ 2.
1145impl<const N: usize> Simple for SUN<N> {}
1146
1147/// SU(N) is semi-simple (implied by simple) for all N ≥ 2.
1148impl<const N: usize> SemiSimple for SUN<N> {}
1149
1150/// su(N) algebra elements are traceless by construction.
1151///
1152/// The representation `SunAlgebra<N>` stores N²-1 coefficients in a
1153/// generalized Gell-Mann basis. All generators are traceless by definition.
1154impl<const N: usize> TracelessByConstruction for SunAlgebra<N> {}
1155
1156/// su(N) algebra elements are anti-Hermitian by construction.
1157///
1158/// The representation uses i·λⱼ where λⱼ are Hermitian generators.
1159impl<const N: usize> AntiHermitianByConstruction for SunAlgebra<N> {}
1160
1161// ============================================================================
1162// Type Aliases
1163// ============================================================================
1164
1165/// Type alias for SU(2) via generic implementation
1166pub type SU2Generic = SUN<2>;
1167/// Type alias for SU(3) via generic implementation
1168pub type SU3Generic = SUN<3>;
1169/// Type alias for SU(4) - Pati-Salam model
1170pub type SU4 = SUN<4>;
1171/// Type alias for SU(5) - Georgi-Glashow GUT
1172pub type SU5 = SUN<5>;
1173
1174#[cfg(test)]
1175mod tests {
1176    use super::*;
1177    use approx::assert_relative_eq;
1178    // ========================================================================
1179    // Cross-implementation consistency: SUN<2> vs SU2, SUN<3> vs SU3
1180    // ========================================================================
1181
1182    /// SPEC: Same coefficients in Su2Algebra and SunAlgebra<2> must produce
1183    /// the same matrix, the same bracket, and the same exp.
1184    ///
1185    /// This is the regression guard for basis normalization consistency.
1186    /// Convention: tr(Tₐ†Tᵦ) = ½δₐᵦ for all implementations.
1187    #[test]
1188    fn test_sun2_su2_basis_matrices_agree() {
1189        // Basis element matrices must be identical
1190        for k in 0..3 {
1191            let m_sun = SunAlgebra::<2>::basis_element(k).to_matrix();
1192            // Build SU2 basis matrix manually: Tₖ = iσₖ/2
1193            let i = Complex64::new(0.0, 1.0);
1194            let sigma: Array2<Complex64> = match k {
1195                0 => Array2::from_shape_vec(
1196                    (2, 2),
1197                    vec![
1198                        Complex64::new(0.0, 0.0),
1199                        Complex64::new(1.0, 0.0),
1200                        Complex64::new(1.0, 0.0),
1201                        Complex64::new(0.0, 0.0),
1202                    ],
1203                )
1204                .unwrap(),
1205                1 => Array2::from_shape_vec(
1206                    (2, 2),
1207                    vec![
1208                        Complex64::new(0.0, 0.0),
1209                        Complex64::new(0.0, -1.0),
1210                        Complex64::new(0.0, 1.0),
1211                        Complex64::new(0.0, 0.0),
1212                    ],
1213                )
1214                .unwrap(),
1215                2 => Array2::from_shape_vec(
1216                    (2, 2),
1217                    vec![
1218                        Complex64::new(1.0, 0.0),
1219                        Complex64::new(0.0, 0.0),
1220                        Complex64::new(0.0, 0.0),
1221                        Complex64::new(-1.0, 0.0),
1222                    ],
1223                )
1224                .unwrap(),
1225                _ => unreachable!(),
1226            };
1227            let expected = sigma.mapv(|z| i * z * 0.5);
1228
1229            for r in 0..2 {
1230                for c in 0..2 {
1231                    assert!(
1232                        (m_sun[(r, c)] - expected[(r, c)]).norm() < 1e-10,
1233                        "SUN<2> basis {} at ({},{}): got {}, want {}",
1234                        k,
1235                        r,
1236                        c,
1237                        m_sun[(r, c)],
1238                        expected[(r, c)]
1239                    );
1240                }
1241            }
1242        }
1243    }
1244
1245    #[test]
1246    fn test_sun2_su2_brackets_agree() {
1247        use crate::Su2Algebra;
1248
1249        for i in 0..3 {
1250            for j in 0..3 {
1251                let bracket_su2 = Su2Algebra::basis_element(i)
1252                    .bracket(&Su2Algebra::basis_element(j))
1253                    .to_components();
1254                let bracket_sun = SunAlgebra::<2>::basis_element(i)
1255                    .bracket(&SunAlgebra::<2>::basis_element(j))
1256                    .to_components();
1257
1258                for k in 0..3 {
1259                    assert_relative_eq!(bracket_su2[k], bracket_sun[k], epsilon = 1e-10);
1260                }
1261            }
1262        }
1263    }
1264
1265    #[test]
1266    fn test_sun2_su2_exp_agrees() {
1267        use crate::{Su2Algebra, SU2};
1268
1269        let coeffs = [0.3, -0.2, 0.4];
1270        let g_su2 = SU2::exp(&Su2Algebra::new(coeffs));
1271        let g_sun = SUN::<2>::exp(&SunAlgebra::<2>::from_components(&coeffs));
1272
1273        let m_su2 = g_su2.to_matrix();
1274        let m_sun = g_sun.matrix();
1275
1276        for i in 0..2 {
1277            for j in 0..2 {
1278                assert_relative_eq!(m_su2[i][j].re, m_sun[(i, j)].re, epsilon = 1e-10);
1279                assert_relative_eq!(m_su2[i][j].im, m_sun[(i, j)].im, epsilon = 1e-10);
1280            }
1281        }
1282    }
1283
1284    /// SPEC: SU3 and SUN<3> basis matrices must agree up to basis reordering.
1285    ///
1286    /// SU3 uses standard Gell-Mann ordering (λ₁..λ₈ interleaved by type).
1287    /// SUN<3> groups by type: symmetric, antisymmetric, diagonal.
1288    /// Mapping: SU3 index [0,1,2,3,4,5,6,7] → SUN<3> index [0,3,6,1,4,2,5,7].
1289    #[test]
1290    fn test_sun3_su3_basis_matrices_agree() {
1291        use crate::Su3Algebra;
1292
1293        // SU3 Gell-Mann index → SUN<3> generalized Gell-Mann index
1294        let su3_to_sun: [usize; 8] = [0, 3, 6, 1, 4, 2, 5, 7];
1295
1296        for k in 0..8 {
1297            let m_su3 = Su3Algebra::basis_element(k).to_matrix();
1298            let m_sun = SunAlgebra::<3>::basis_element(su3_to_sun[k]).to_matrix();
1299
1300            for r in 0..3 {
1301                for c in 0..3 {
1302                    assert!(
1303                        (m_su3[(r, c)] - m_sun[(r, c)]).norm() < 1e-10,
1304                        "Basis {} (SU3) vs {} (SUN) at ({},{}): SU3={}, SUN<3>={}",
1305                        k,
1306                        su3_to_sun[k],
1307                        r,
1308                        c,
1309                        m_su3[(r, c)],
1310                        m_sun[(r, c)]
1311                    );
1312                }
1313            }
1314        }
1315    }
1316
1317    /// SPEC: Brackets must agree when using corresponding basis elements.
1318    ///
1319    /// We compare at the matrix level to avoid coefficient ordering issues.
1320    #[test]
1321    fn test_sun3_su3_brackets_agree() {
1322        use crate::Su3Algebra;
1323
1324        // SU3 Gell-Mann index → SUN<3> generalized Gell-Mann index
1325        let su3_to_sun: [usize; 8] = [0, 3, 6, 1, 4, 2, 5, 7];
1326
1327        for i in 0..8 {
1328            for j in 0..8 {
1329                let bracket_su3 = Su3Algebra::basis_element(i)
1330                    .bracket(&Su3Algebra::basis_element(j))
1331                    .to_matrix();
1332                let bracket_sun = SunAlgebra::<3>::basis_element(su3_to_sun[i])
1333                    .bracket(&SunAlgebra::<3>::basis_element(su3_to_sun[j]))
1334                    .to_matrix();
1335
1336                for r in 0..3 {
1337                    for c in 0..3 {
1338                        assert!(
1339                            (bracket_su3[(r, c)] - bracket_sun[(r, c)]).norm() < 1e-10,
1340                            "Bracket [e{},e{}] at ({},{}): SU3={}, SUN<3>={}",
1341                            i,
1342                            j,
1343                            r,
1344                            c,
1345                            bracket_su3[(r, c)],
1346                            bracket_sun[(r, c)]
1347                        );
1348                    }
1349                }
1350            }
1351        }
1352    }
1353
1354    /// SPEC: exp of the same matrix must produce the same group element.
1355    ///
1356    /// We construct the algebra element via matrix to avoid ordering issues.
1357    #[test]
1358    fn test_sun3_su3_exp_agrees() {
1359        use crate::{Su3Algebra, SU3};
1360
1361        // Build an su(3) matrix via SU3, then reconstruct in both representations
1362        let su3_coeffs = [0.1, -0.2, 0.15, 0.08, -0.12, 0.05, 0.1, -0.06];
1363        let su3_elem = Su3Algebra::from_components(&su3_coeffs);
1364        let matrix = su3_elem.to_matrix();
1365
1366        // Reconstruct via SUN<3>::from_matrix
1367        let sun_elem = SunAlgebra::<3>::from_matrix(&matrix);
1368
1369        let g_su3 = SU3::exp(&su3_elem);
1370        let g_sun = SUN::<3>::exp(&sun_elem);
1371
1372        let m_su3 = g_su3.matrix();
1373        let m_sun = g_sun.matrix();
1374
1375        for i in 0..3 {
1376            for j in 0..3 {
1377                assert!(
1378                    (m_su3[(i, j)] - m_sun[(i, j)]).norm() < 1e-10,
1379                    "exp disagrees at ({},{}): SU3={}, SUN<3>={}",
1380                    i,
1381                    j,
1382                    m_su3[(i, j)],
1383                    m_sun[(i, j)]
1384                );
1385            }
1386        }
1387    }
1388
1389    /// SPEC: All implementations must satisfy tr(Tₐ†Tᵦ) = ½δₐᵦ.
1390    /// Tests both diagonal (normalization) and off-diagonal (orthogonality).
1391    #[test]
1392    fn test_normalization_half_delta() {
1393        for n in [2_usize, 3, 4, 5] {
1394            let dim = n * n - 1;
1395            for a in 0..dim {
1396                let ma = match n {
1397                    2 => SunAlgebra::<2>::basis_element(a).to_matrix(),
1398                    3 => SunAlgebra::<3>::basis_element(a).to_matrix(),
1399                    4 => SunAlgebra::<4>::basis_element(a).to_matrix(),
1400                    5 => SunAlgebra::<5>::basis_element(a).to_matrix(),
1401                    _ => unreachable!(),
1402                };
1403                let ma_dag = ma.t().mapv(|z| z.conj());
1404
1405                for b in a..dim {
1406                    let mb = match n {
1407                        2 => SunAlgebra::<2>::basis_element(b).to_matrix(),
1408                        3 => SunAlgebra::<3>::basis_element(b).to_matrix(),
1409                        4 => SunAlgebra::<4>::basis_element(b).to_matrix(),
1410                        5 => SunAlgebra::<5>::basis_element(b).to_matrix(),
1411                        _ => unreachable!(),
1412                    };
1413
1414                    let prod = ma_dag.dot(&mb);
1415                    let mut tr = 0.0;
1416                    for i in 0..n {
1417                        tr += prod[(i, i)].re;
1418                    }
1419
1420                    let expected = if a == b { 0.5 } else { 0.0 };
1421                    assert!(
1422                        (tr - expected).abs() < 1e-10,
1423                        "SU({}): tr(T{}†T{}) = {:.4}, want {}",
1424                        n,
1425                        a,
1426                        b,
1427                        tr,
1428                        expected
1429                    );
1430                }
1431            }
1432        }
1433    }
1434
1435    // ========================================================================
1436    // Normalization consistency: bracket + exp + adjoint must be coherent
1437    // ========================================================================
1438
1439    /// Verify BCH(X,Y) ≈ log(exp(X)·exp(Y)) for SUN<2>, confirming
1440    /// the bracket/exp coupling is correct under the /2 convention.
1441    #[test]
1442    fn test_bch_exp_log_coherence_sun2() {
1443        use crate::bch::bch_fifth_order;
1444
1445        let x = SunAlgebra::<2>::from_components(&[0.05, -0.03, 0.04]);
1446        let y = SunAlgebra::<2>::from_components(&[-0.02, 0.04, -0.03]);
1447
1448        let direct = SUN::<2>::exp(&x).compose(&SUN::<2>::exp(&y));
1449        let bch_z = bch_fifth_order(&x, &y);
1450        let via_bch = SUN::<2>::exp(&bch_z);
1451
1452        let distance = direct.compose(&via_bch.inverse()).distance_to_identity();
1453        assert!(
1454            distance < 1e-8,
1455            "SUN<2> BCH vs exp·log distance: {:.2e}",
1456            distance
1457        );
1458    }
1459
1460    /// Adjoint preserves bracket on SU3: Ad_g([X,Y]) = [Ad_g(X), Ad_g(Y)].
1461    /// Uses non-random, non-axis-aligned elements for deterministic coverage.
1462    #[test]
1463    fn test_adjoint_preserves_bracket_su3() {
1464        use crate::{Su3Algebra, SU3};
1465
1466        let x = Su3Algebra::from_components(&[0.3, -0.2, 0.1, 0.15, -0.1, 0.25, -0.15, 0.05]);
1467        let y = Su3Algebra::from_components(&[-0.1, 0.25, -0.15, 0.2, 0.05, -0.1, 0.3, -0.2]);
1468        let g = SU3::exp(&Su3Algebra::from_components(&[0.4, -0.3, 0.2, 0.1, -0.2, 0.15, -0.1, 0.3]));
1469
1470        let lhs = g.adjoint_action(&x.bracket(&y));
1471        let rhs = g.adjoint_action(&x).bracket(&g.adjoint_action(&y));
1472
1473        let diff_matrix = lhs.to_matrix() - rhs.to_matrix();
1474        let mut frobenius = 0.0;
1475        for r in 0..3 {
1476            for c in 0..3 {
1477                frobenius += diff_matrix[(r, c)].norm_sqr();
1478            }
1479        }
1480        assert!(
1481            frobenius.sqrt() < 1e-10,
1482            "SU3 Ad_g([X,Y]) ≠ [Ad_g(X), Ad_g(Y)]: ||diff|| = {:.2e}",
1483            frobenius.sqrt()
1484        );
1485    }
1486
1487    /// SU3 bracket (structure constants) and SUN<3> bracket (matrix commutator)
1488    /// must produce the same abstract element, verified through from_matrix roundtrip.
1489    #[test]
1490    fn test_su3_sun3_bracket_coefficient_roundtrip() {
1491        use crate::Su3Algebra;
1492
1493        let x_su3 = Su3Algebra::from_components(&[0.3, -0.2, 0.1, 0.15, -0.1, 0.25, -0.15, 0.05]);
1494        let y_su3 = Su3Algebra::from_components(&[-0.1, 0.25, -0.15, 0.2, 0.05, -0.1, 0.3, -0.2]);
1495
1496        // Bracket via SU3 structure constants
1497        let bracket_su3 = x_su3.bracket(&y_su3);
1498        let m_su3 = bracket_su3.to_matrix();
1499
1500        // Same elements via SUN<3>, bracket via matrix commutator
1501        let x_sun = SunAlgebra::<3>::from_matrix(&x_su3.to_matrix());
1502        let y_sun = SunAlgebra::<3>::from_matrix(&y_su3.to_matrix());
1503        let bracket_sun = x_sun.bracket(&y_sun);
1504        let m_sun = bracket_sun.to_matrix();
1505
1506        for r in 0..3 {
1507            for c in 0..3 {
1508                assert!(
1509                    (m_su3[(r, c)] - m_sun[(r, c)]).norm() < 1e-12,
1510                    "Bracket matrix disagrees at ({},{}): SU3={}, SUN<3>={}",
1511                    r, c, m_su3[(r, c)], m_sun[(r, c)]
1512                );
1513            }
1514        }
1515
1516        // Also verify from_matrix roundtrip: SU3 → matrix → SUN<3> → matrix
1517        let roundtrip = SunAlgebra::<3>::from_matrix(&m_su3).to_matrix();
1518        for r in 0..3 {
1519            for c in 0..3 {
1520                assert!(
1521                    (m_su3[(r, c)] - roundtrip[(r, c)]).norm() < 1e-12,
1522                    "from_matrix roundtrip failed at ({},{})",
1523                    r, c
1524                );
1525            }
1526        }
1527    }
1528
1529    // ========================================================================
1530    // Group axioms for N=4,5 (untested dimensions)
1531    // ========================================================================
1532
1533    #[test]
1534    fn test_su4_group_axioms() {
1535        let x = SunAlgebra::<4>::from_components(&[
1536            0.1, -0.2, 0.15, 0.08, -0.12, 0.05, 0.1, -0.06, 0.09, -0.11, 0.07, 0.03, -0.08, 0.04,
1537            0.13,
1538        ]);
1539        let y = SunAlgebra::<4>::from_components(&[
1540            -0.05, 0.1, -0.08, 0.12, 0.06, -0.15, 0.03, 0.09, -0.07, 0.11, -0.04, 0.14, 0.02, -0.1,
1541            0.08,
1542        ]);
1543
1544        let g = SUN::<4>::exp(&x);
1545        let h = SUN::<4>::exp(&y);
1546        let e = SUN::<4>::identity();
1547
1548        // Identity
1549        assert_relative_eq!(
1550            g.compose(&e).distance_to_identity(),
1551            g.distance_to_identity(),
1552            epsilon = 1e-10
1553        );
1554
1555        // Inverse
1556        assert_relative_eq!(
1557            g.compose(&g.inverse()).distance_to_identity(),
1558            0.0,
1559            epsilon = 1e-9
1560        );
1561
1562        // Associativity
1563        let k = SUN::<4>::exp(&x.scale(0.7));
1564        let lhs = g.compose(&h).compose(&k);
1565        let rhs = g.compose(&h.compose(&k));
1566        assert_relative_eq!(
1567            lhs.compose(&rhs.inverse()).distance_to_identity(),
1568            0.0,
1569            epsilon = 1e-9
1570        );
1571
1572        // Unitarity preservation
1573        assert!(g.verify_unitarity(1e-10));
1574        assert!(g.compose(&h).verify_unitarity(1e-10));
1575    }
1576
1577    #[test]
1578    fn test_su5_group_axioms() {
1579        let x = SunAlgebra::<5>::from_components(
1580            &(0..24)
1581                .map(|i| 0.05 * (i as f64 - 12.0).sin())
1582                .collect::<Vec<_>>(),
1583        );
1584        let g = SUN::<5>::exp(&x);
1585
1586        // Identity
1587        assert_relative_eq!(
1588            g.compose(&SUN::<5>::identity()).distance_to_identity(),
1589            g.distance_to_identity(),
1590            epsilon = 1e-10
1591        );
1592
1593        // Inverse
1594        assert_relative_eq!(
1595            g.compose(&g.inverse()).distance_to_identity(),
1596            0.0,
1597            epsilon = 1e-9
1598        );
1599
1600        // Unitarity
1601        assert!(g.verify_unitarity(1e-10));
1602    }
1603
1604    // ========================================================================
1605    // Jacobi identity for N=4 (tested only for N=3 previously)
1606    // ========================================================================
1607
1608    #[test]
1609    fn test_su4_jacobi_identity() {
1610        // Test with several triples of basis elements
1611        let triples = [(0, 3, 7), (1, 5, 10), (2, 8, 14), (4, 9, 12)];
1612        for (i, j, k) in triples {
1613            let x = SunAlgebra::<4>::basis_element(i);
1614            let y = SunAlgebra::<4>::basis_element(j);
1615            let z = SunAlgebra::<4>::basis_element(k);
1616
1617            let t1 = x.bracket(&y.bracket(&z));
1618            let t2 = y.bracket(&z.bracket(&x));
1619            let t3 = z.bracket(&x.bracket(&y));
1620            let sum = t1.add(&t2).add(&t3);
1621
1622            assert!(
1623                sum.norm() < 1e-10,
1624                "Jacobi violated for SU(4) basis ({},{},{}): ||sum|| = {:.2e}",
1625                i,
1626                j,
1627                k,
1628                sum.norm()
1629            );
1630        }
1631    }
1632
1633    // ========================================================================
1634    // Adjoint representation axioms for SU(N)
1635    // ========================================================================
1636
1637    #[test]
1638    fn test_sun_adjoint_identity_action() {
1639        // Ad_e(X) = X for all X
1640        let e = SUN::<3>::identity();
1641        for i in 0..8 {
1642            let x = SunAlgebra::<3>::basis_element(i);
1643            let ad_x = e.adjoint_action(&x);
1644            for k in 0..8 {
1645                assert_relative_eq!(
1646                    x.to_components()[k],
1647                    ad_x.to_components()[k],
1648                    epsilon = 1e-10
1649                );
1650            }
1651        }
1652    }
1653
1654    #[test]
1655    fn test_sun_adjoint_homomorphism() {
1656        // Ad_{g·h}(X) = Ad_g(Ad_h(X))
1657        let g = SUN::<3>::exp(&SunAlgebra::<3>::basis_element(0).scale(0.5));
1658        let h = SUN::<3>::exp(&SunAlgebra::<3>::basis_element(3).scale(0.3));
1659        let x = SunAlgebra::<3>::from_components(&[0.1, -0.2, 0.15, 0.08, -0.12, 0.05, 0.1, -0.06]);
1660
1661        let gh = g.compose(&h);
1662        let lhs = gh.adjoint_action(&x);
1663        let rhs = g.adjoint_action(&h.adjoint_action(&x));
1664
1665        for k in 0..8 {
1666            assert_relative_eq!(
1667                lhs.to_components()[k],
1668                rhs.to_components()[k],
1669                epsilon = 1e-9
1670            );
1671        }
1672    }
1673
1674    #[test]
1675    fn test_sun_adjoint_preserves_bracket() {
1676        // Ad_g([X,Y]) = [Ad_g(X), Ad_g(Y)]
1677        let g = SUN::<3>::exp(&SunAlgebra::<3>::basis_element(2).scale(0.8));
1678        let x = SunAlgebra::<3>::basis_element(0);
1679        let y = SunAlgebra::<3>::basis_element(4);
1680
1681        let lhs = g.adjoint_action(&x.bracket(&y));
1682        let rhs = g.adjoint_action(&x).bracket(&g.adjoint_action(&y));
1683
1684        for k in 0..8 {
1685            assert_relative_eq!(
1686                lhs.to_components()[k],
1687                rhs.to_components()[k],
1688                epsilon = 1e-9
1689            );
1690        }
1691    }
1692
1693    #[test]
1694    fn test_sun4_adjoint_preserves_norm() {
1695        // ||Ad_g(X)|| = ||X|| for compact groups
1696        let x = SunAlgebra::<4>::from_components(&[
1697            0.1, -0.2, 0.15, 0.08, -0.12, 0.05, 0.1, -0.06, 0.09, -0.11, 0.07, 0.03, -0.08, 0.04,
1698            0.13,
1699        ]);
1700        let g = SUN::<4>::exp(&SunAlgebra::<4>::basis_element(5).scale(1.2));
1701        let ad_x = g.adjoint_action(&x);
1702        assert_relative_eq!(x.norm(), ad_x.norm(), epsilon = 1e-9);
1703    }
1704
1705    // ========================================================================
1706    // exp/log roundtrip for N=4 (untested dimension)
1707    // ========================================================================
1708
1709    #[test]
1710    fn test_su4_exp_log_roundtrip() {
1711        let x = SunAlgebra::<4>::from_components(&[
1712            0.1, -0.05, 0.08, 0.03, -0.06, 0.02, 0.04, -0.07, 0.09, -0.03, 0.01, 0.05, -0.02, 0.06,
1713            0.04,
1714        ]);
1715        let g = SUN::<4>::exp(&x);
1716        assert!(g.verify_unitarity(1e-10));
1717
1718        let x_back = SUN::<4>::log(&g).expect("log should succeed near identity");
1719        let diff: f64 = x
1720            .coefficients()
1721            .iter()
1722            .zip(x_back.coefficients().iter())
1723            .map(|(a, b)| (a - b).powi(2))
1724            .sum::<f64>()
1725            .sqrt();
1726
1727        assert!(diff < 1e-8, "SU(4) exp/log roundtrip error: {:.2e}", diff);
1728    }
1729
1730    #[test]
1731    fn test_sun_algebra_dimensions() {
1732        assert_eq!(SunAlgebra::<2>::DIM, 3); // SU(2)
1733        assert_eq!(SunAlgebra::<3>::DIM, 8); // SU(3)
1734        assert_eq!(SunAlgebra::<4>::DIM, 15); // SU(4)
1735        assert_eq!(SunAlgebra::<5>::DIM, 24); // SU(5)
1736    }
1737
1738    #[test]
1739    fn test_sun_algebra_zero() {
1740        let zero = SunAlgebra::<3>::zero();
1741        assert_eq!(zero.coefficients.len(), 8);
1742        assert!(zero.coefficients.iter().all(|&x| x == 0.0));
1743    }
1744
1745    #[test]
1746    fn test_sun_algebra_add_scale() {
1747        let x = SunAlgebra::<2>::basis_element(0);
1748        let y = SunAlgebra::<2>::basis_element(1);
1749
1750        let sum = &x + &y;
1751        assert_eq!(sum.coefficients, vec![1.0, 1.0, 0.0]);
1752
1753        let scaled = x.scale(2.5);
1754        assert_eq!(scaled.coefficients, vec![2.5, 0.0, 0.0]);
1755    }
1756
1757    #[test]
1758    fn test_sun_identity() {
1759        let id = SUN::<3>::identity();
1760        assert!(id.verify_unitarity(1e-10));
1761        assert_relative_eq!(id.distance_to_identity(), 0.0, epsilon = 1e-10);
1762    }
1763
1764    #[test]
1765    fn test_sun_exponential_preserves_unitarity() {
1766        // Random algebra element
1767        let algebra =
1768            SunAlgebra::<3>::from_components(&[0.5, -0.3, 0.8, 0.2, -0.6, 0.4, 0.1, -0.2]);
1769        let g = SUN::<3>::exp(&algebra);
1770
1771        // Verify U†U = I
1772        assert!(
1773            g.verify_unitarity(1e-10),
1774            "Exponential should preserve unitarity"
1775        );
1776    }
1777
1778    #[test]
1779    fn test_sun_exp_identity() {
1780        let zero = SunAlgebra::<4>::zero();
1781        let g = SUN::<4>::exp(&zero);
1782        assert_relative_eq!(g.distance_to_identity(), 0.0, epsilon = 1e-10);
1783    }
1784
1785    #[test]
1786    fn test_sun_group_composition() {
1787        let g1 = SUN::<2>::exp(&SunAlgebra::<2>::basis_element(0).scale(0.5));
1788        let g2 = SUN::<2>::exp(&SunAlgebra::<2>::basis_element(1).scale(0.3));
1789
1790        let product = g1.compose(&g2);
1791
1792        assert!(product.verify_unitarity(1e-10));
1793    }
1794
1795    #[test]
1796    fn test_sun_inverse() {
1797        let algebra =
1798            SunAlgebra::<3>::from_components(&[0.2, 0.3, -0.1, 0.5, -0.2, 0.1, 0.4, -0.3]);
1799        let g = SUN::<3>::exp(&algebra);
1800        let g_inv = g.inverse();
1801
1802        let product = g.compose(&g_inv);
1803
1804        assert_relative_eq!(product.distance_to_identity(), 0.0, epsilon = 1e-9);
1805    }
1806
1807    #[test]
1808    fn test_sun_adjoint_action_preserves_norm() {
1809        let g = SUN::<3>::exp(&SunAlgebra::<3>::basis_element(0).scale(1.2));
1810        let x = SunAlgebra::<3>::basis_element(2).scale(0.5);
1811
1812        let ad_x = g.adjoint_action(&x);
1813
1814        // Adjoint action preserves norm for compact groups
1815        assert_relative_eq!(x.norm(), ad_x.norm(), epsilon = 1e-9);
1816    }
1817
1818    #[test]
1819    fn test_sun_exp_log_roundtrip() {
1820        // SU(3): exp then log should recover the original algebra element
1821        let x = SunAlgebra::<3>::from_components(&[0.1, -0.2, 0.15, 0.08, -0.12, 0.05, 0.1, -0.06]);
1822        let g = SUN::<3>::exp(&x);
1823        assert!(g.verify_unitarity(1e-10));
1824
1825        let x_back = SUN::<3>::log(&g).expect("log should succeed near identity");
1826        let diff_norm: f64 = x
1827            .coefficients
1828            .iter()
1829            .zip(x_back.coefficients.iter())
1830            .map(|(a, b)| (a - b).powi(2))
1831            .sum::<f64>()
1832            .sqrt();
1833
1834        assert!(
1835            diff_norm < 1e-8,
1836            "SU(3) exp/log roundtrip error: {:.2e}",
1837            diff_norm
1838        );
1839    }
1840
1841    #[test]
1842    fn test_sun_exp_log_roundtrip_su2() {
1843        // SU(2) via generic SUN: smaller algebra
1844        let x = SunAlgebra::<2>::from_components(&[0.3, -0.2, 0.4]);
1845        let g = SUN::<2>::exp(&x);
1846        assert!(g.verify_unitarity(1e-10));
1847
1848        let x_back = SUN::<2>::log(&g).expect("log should succeed");
1849        let diff_norm: f64 = x
1850            .coefficients
1851            .iter()
1852            .zip(x_back.coefficients.iter())
1853            .map(|(a, b)| (a - b).powi(2))
1854            .sum::<f64>()
1855            .sqrt();
1856
1857        assert!(
1858            diff_norm < 1e-8,
1859            "SU(2) exp/log roundtrip error: {:.2e}",
1860            diff_norm
1861        );
1862    }
1863
1864    #[test]
1865    fn test_sun_log_exp_roundtrip() {
1866        // Start from group element, log, then exp back
1867        let x = SunAlgebra::<3>::from_components(&[0.2, 0.3, -0.1, 0.5, -0.2, 0.1, 0.4, -0.3]);
1868        let g = SUN::<3>::exp(&x);
1869
1870        let log_g = SUN::<3>::log(&g).expect("log should succeed");
1871        let g_back = SUN::<3>::exp(&log_g);
1872
1873        assert_relative_eq!(
1874            g.distance_to_identity(),
1875            g_back.distance_to_identity(),
1876            epsilon = 1e-8
1877        );
1878
1879        // Check that g and g_back are close
1880        let product = g.compose(&g_back.inverse());
1881        assert_relative_eq!(product.distance_to_identity(), 0.0, epsilon = 1e-8);
1882    }
1883
1884    #[test]
1885    fn test_sun_jacobi_identity() {
1886        let x = SunAlgebra::<3>::basis_element(0);
1887        let y = SunAlgebra::<3>::basis_element(1);
1888        let z = SunAlgebra::<3>::basis_element(2);
1889
1890        // [X, [Y, Z]] + [Y, [Z, X]] + [Z, [X, Y]] = 0
1891        let t1 = x.bracket(&y.bracket(&z));
1892        let t2 = y.bracket(&z.bracket(&x));
1893        let t3 = z.bracket(&x.bracket(&y));
1894        let sum = t1.add(&t2).add(&t3);
1895
1896        assert!(
1897            sum.norm() < 1e-10,
1898            "Jacobi identity violated for SU(3): ||sum|| = {:.2e}",
1899            sum.norm()
1900        );
1901    }
1902
1903    #[test]
1904    fn test_sun_bracket_antisymmetry() {
1905        let x = SunAlgebra::<4>::basis_element(0);
1906        let y = SunAlgebra::<4>::basis_element(3);
1907
1908        let xy = x.bracket(&y);
1909        let yx = y.bracket(&x);
1910
1911        // [X,Y] = -[Y,X]
1912        for i in 0..SunAlgebra::<4>::DIM {
1913            assert_relative_eq!(xy.coefficients[i], -yx.coefficients[i], epsilon = 1e-10);
1914        }
1915    }
1916
1917    #[test]
1918    fn test_sun_bracket_bilinearity() {
1919        let x = SunAlgebra::<3>::basis_element(0);
1920        let y = SunAlgebra::<3>::basis_element(3);
1921        let z = SunAlgebra::<3>::basis_element(5);
1922        let alpha = 2.5;
1923
1924        // Left linearity: [αX + Y, Z] = α[X, Z] + [Y, Z]
1925        let lhs = x.scale(alpha).add(&y).bracket(&z);
1926        let rhs = x.bracket(&z).scale(alpha).add(&y.bracket(&z));
1927        for i in 0..SunAlgebra::<3>::DIM {
1928            assert!(
1929                (lhs.coefficients[i] - rhs.coefficients[i]).abs() < 1e-14,
1930                "Left linearity failed at component {}: {} vs {}",
1931                i,
1932                lhs.coefficients[i],
1933                rhs.coefficients[i]
1934            );
1935        }
1936
1937        // Right linearity: [Z, αX + Y] = α[Z, X] + [Z, Y]
1938        let lhs = z.bracket(&x.scale(alpha).add(&y));
1939        let rhs = z.bracket(&x).scale(alpha).add(&z.bracket(&y));
1940        for i in 0..SunAlgebra::<3>::DIM {
1941            assert!(
1942                (lhs.coefficients[i] - rhs.coefficients[i]).abs() < 1e-14,
1943                "Right linearity failed at component {}: {} vs {}",
1944                i,
1945                lhs.coefficients[i],
1946                rhs.coefficients[i]
1947            );
1948        }
1949    }
1950}