Skip to main content

lie_groups/
su3.rs

1//! Lie group SU(3) - Special unitary 3×3 group
2//!
3//! SU(3) is the group of 3×3 complex unitary matrices with determinant 1.
4//! It is the gauge group of quantum chromodynamics (QCD).
5//!
6//! # Mathematical Structure
7//!
8//! ```text
9//! SU(3) = { U ∈ ℂ³ˣ³ | U† U = I, det(U) = 1 }
10//! ```
11//!
12//! # Lie Algebra
13//!
14//! The Lie algebra su(3) consists of 3×3 traceless anti-Hermitian matrices:
15//! ```text
16//! su(3) = { X ∈ ℂ³ˣ³ | X† = -X, Tr(X) = 0 }
17//! ```
18//!
19//! This is 8-dimensional, with basis given by the Gell-Mann matrices.
20//!
21//! # Gell-Mann Matrices
22//!
23//! The 8 generators λ₁, ..., λ₈ are:
24//! ```text
25//! λ₁ = [[0,1,0],[1,0,0],[0,0,0]]           (like Pauli σₓ in 1-2 block)
26//! λ₂ = [[0,-i,0],[i,0,0],[0,0,0]]          (like Pauli σᵧ in 1-2 block)
27//! λ₃ = [[1,0,0],[0,-1,0],[0,0,0]]          (like Pauli σᵤ in 1-2 block)
28//! λ₄ = [[0,0,1],[0,0,0],[1,0,0]]           (1-3 mixing)
29//! λ₅ = [[0,0,-i],[0,0,0],[i,0,0]]          (1-3 mixing)
30//! λ₆ = [[0,0,0],[0,0,1],[0,1,0]]           (2-3 mixing)
31//! λ₇ = [[0,0,0],[0,0,-i],[0,i,0]]          (2-3 mixing)
32//! λ₈ = (1/√3)[[1,0,0],[0,1,0],[0,0,-2]]    (diagonal, traceless)
33//! ```
34//!
35//! # Color Charge
36//!
37//! In QCD, SU(3) acts on the 3-dimensional color space (red, green, blue).
38//! Quarks transform in the fundamental representation, gluons in the adjoint.
39//!
40//! # Scaling-and-Squaring Algorithm
41//!
42//! The matrix exponential uses scaling-and-squaring (Higham, "Functions of Matrices"):
43//!
44//! ```text
45//! exp(X) = [exp(X/2^k)]^{2^k}
46//! ```
47//!
48//! ## Threshold Choice: ||X/2^k|| ≤ 0.5
49//!
50//! We choose k such that the scaled matrix has Frobenius norm ≤ 0.5. This threshold
51//! balances accuracy and efficiency:
52//!
53//! **Convergence Analysis:**
54//! The Taylor series exp(Y) = I + Y + Y²/2! + ... has remainder bounded by:
55//! ```text
56//! ||exp(Y) - Σₙ Yⁿ/n!|| ≤ ||Y||^{N+1} / (N+1)! × exp(||Y||)
57//! ```
58//!
59//! For ||Y|| ≤ 0.5 and N = 15 terms:
60//! - Truncation error: 0.5^16 / 16! × e^0.5 ≈ 4 × 10^{-18}
61//! - Well below IEEE 754 f64 machine epsilon (2.2 × 10^{-16})
62//!
63//! **Why not smaller?** Using ||Y|| ≤ 0.1 requires more squarings (k larger),
64//! and each squaring accumulates rounding error. The choice 0.5 minimizes total error.
65//!
66//! **Why not larger?** For ||Y|| > 1, the series converges slowly and requires
67//! many more terms, increasing both computation and accumulated rounding error.
68//!
69//! **Reference:** Al-Mohy & Higham, "A New Scaling and Squaring Algorithm for the
70//! Matrix Exponential" (2010), SIAM J. Matrix Anal. Appl.
71
72use crate::traits::{AntiHermitianByConstruction, LieAlgebra, LieGroup, TracelessByConstruction};
73use ndarray::Array2;
74use num_complex::Complex64;
75use std::fmt;
76use std::ops::{Add, Mul, MulAssign, Neg, Sub};
77
78/// Lie algebra su(3) - 8-dimensional space of traceless anti-Hermitian 3×3 matrices
79///
80/// Elements are represented by 8 real coefficients [a₁, a₂, ..., a₈] corresponding
81/// to the linear combination:
82/// ```text
83/// X = i·∑ⱼ aⱼ·λⱼ
84/// ```
85/// where λⱼ are the Gell-Mann matrices (j = 1..8).
86///
87/// # Basis Elements
88///
89/// The 8 Gell-Mann matrices form a basis for su(3). They satisfy:
90/// - Hermitian: λⱼ† = λⱼ
91/// - Traceless: Tr(λⱼ) = 0
92/// - Normalized: Tr(λⱼ λₖ) = 2δⱼₖ
93///
94/// # Examples
95///
96/// ```
97/// use lie_groups::su3::Su3Algebra;
98/// use lie_groups::traits::LieAlgebra;
99///
100/// // First basis element (λ₁)
101/// let e1 = Su3Algebra::basis_element(0);
102/// assert_eq!(e1.0, [1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]);
103/// ```
104#[derive(Clone, Copy, Debug, PartialEq)]
105pub struct Su3Algebra(pub [f64; 8]);
106
107impl Add for Su3Algebra {
108    type Output = Self;
109    fn add(self, rhs: Self) -> Self {
110        let mut r = [0.0; 8];
111        for i in 0..8 {
112            r[i] = self.0[i] + rhs.0[i];
113        }
114        Self(r)
115    }
116}
117
118impl Add<&Su3Algebra> for Su3Algebra {
119    type Output = Su3Algebra;
120    fn add(self, rhs: &Su3Algebra) -> Su3Algebra {
121        self + *rhs
122    }
123}
124
125impl Add<Su3Algebra> for &Su3Algebra {
126    type Output = Su3Algebra;
127    fn add(self, rhs: Su3Algebra) -> Su3Algebra {
128        *self + rhs
129    }
130}
131
132impl Add<&Su3Algebra> for &Su3Algebra {
133    type Output = Su3Algebra;
134    fn add(self, rhs: &Su3Algebra) -> Su3Algebra {
135        *self + *rhs
136    }
137}
138
139impl Sub for Su3Algebra {
140    type Output = Self;
141    fn sub(self, rhs: Self) -> Self {
142        let mut r = [0.0; 8];
143        for i in 0..8 {
144            r[i] = self.0[i] - rhs.0[i];
145        }
146        Self(r)
147    }
148}
149
150impl Neg for Su3Algebra {
151    type Output = Self;
152    fn neg(self) -> Self {
153        let mut r = [0.0; 8];
154        for i in 0..8 {
155            r[i] = -self.0[i];
156        }
157        Self(r)
158    }
159}
160
161impl Mul<f64> for Su3Algebra {
162    type Output = Self;
163    fn mul(self, scalar: f64) -> Self {
164        let mut r = [0.0; 8];
165        for i in 0..8 {
166            r[i] = self.0[i] * scalar;
167        }
168        Self(r)
169    }
170}
171
172impl Mul<Su3Algebra> for f64 {
173    type Output = Su3Algebra;
174    fn mul(self, rhs: Su3Algebra) -> Su3Algebra {
175        rhs * self
176    }
177}
178
179impl LieAlgebra for Su3Algebra {
180    const DIM: usize = 8;
181
182    fn zero() -> Self {
183        Self([0.0; 8])
184    }
185
186    fn add(&self, other: &Self) -> Self {
187        let mut result = [0.0; 8];
188        for i in 0..8 {
189            result[i] = self.0[i] + other.0[i];
190        }
191        Self(result)
192    }
193
194    fn scale(&self, scalar: f64) -> Self {
195        let mut result = [0.0; 8];
196        for i in 0..8 {
197            result[i] = self.0[i] * scalar;
198        }
199        Self(result)
200    }
201
202    fn norm(&self) -> f64 {
203        self.0.iter().map(|x| x.powi(2)).sum::<f64>().sqrt()
204    }
205
206    fn basis_element(i: usize) -> Self {
207        assert!(i < 8, "SU(3) algebra is 8-dimensional");
208        let mut coeffs = [0.0; 8];
209        coeffs[i] = 1.0;
210        Self(coeffs)
211    }
212
213    fn from_components(components: &[f64]) -> Self {
214        assert_eq!(components.len(), 8, "su(3) has dimension 8");
215        let mut coeffs = [0.0; 8];
216        coeffs.copy_from_slice(components);
217        Self(coeffs)
218    }
219
220    fn to_components(&self) -> Vec<f64> {
221        self.0.to_vec()
222    }
223
224    /// Lie bracket for su(3): [X, Y] = XY - YX
225    ///
226    /// Computed using structure constants fᵢⱼₖ for O(1) performance.
227    ///
228    /// # Mathematical Formula
229    ///
230    /// For X = i·∑ᵢ aᵢ·λᵢ and Y = i·∑ⱼ bⱼ·λⱼ:
231    /// ```text
232    /// [X, Y] = -2i·∑ₖ (∑ᵢⱼ aᵢ·bⱼ·fᵢⱼₖ)·λₖ
233    /// ```
234    ///
235    /// # Performance
236    ///
237    /// Uses pre-computed table of non-zero structure constants.
238    /// - Old implementation: O(512) iterations with conditional checks
239    /// - New implementation: O(54) direct lookups
240    /// - Speedup: ~10× fewer operations
241    ///
242    /// # Complexity
243    ///
244    /// O(1) - constant time (54 multiply-adds)
245    ///
246    /// # Properties
247    ///
248    /// - Antisymmetric: [X, Y] = -[Y, X]
249    /// - Jacobi identity: [X, [Y, Z]] + [Y, [Z, X]] + [Z, [X, Y]] = 0
250    fn bracket(&self, other: &Self) -> Self {
251        // Pre-computed non-zero structure constants f_ijk
252        // SU(3) has 9 fundamental non-zero f values, each with 6 permutations = 54 total
253        // Format: (i, j, k, f_ijk)
254        const SQRT3_HALF: f64 = 0.866_025_403_784_438_6;
255
256        // All 54 non-zero structure constant entries
257        // Grouped by fundamental: 3 even perms (+f) + 3 odd perms (-f)
258        #[rustfmt::skip]
259        const STRUCTURE_CONSTANTS: [(usize, usize, usize, f64); 54] = [
260            // f_012 = 1
261            (0, 1, 2, 1.0), (1, 2, 0, 1.0), (2, 0, 1, 1.0),
262            (1, 0, 2, -1.0), (0, 2, 1, -1.0), (2, 1, 0, -1.0),
263            // f_036 = 0.5
264            (0, 3, 6, 0.5), (3, 6, 0, 0.5), (6, 0, 3, 0.5),
265            (3, 0, 6, -0.5), (0, 6, 3, -0.5), (6, 3, 0, -0.5),
266            // f_045 = -0.5
267            (0, 4, 5, -0.5), (4, 5, 0, -0.5), (5, 0, 4, -0.5),
268            (4, 0, 5, 0.5), (0, 5, 4, 0.5), (5, 4, 0, 0.5),
269            // f_135 = 0.5
270            (1, 3, 5, 0.5), (3, 5, 1, 0.5), (5, 1, 3, 0.5),
271            (3, 1, 5, -0.5), (1, 5, 3, -0.5), (5, 3, 1, -0.5),
272            // f_146 = 0.5
273            (1, 4, 6, 0.5), (4, 6, 1, 0.5), (6, 1, 4, 0.5),
274            (4, 1, 6, -0.5), (1, 6, 4, -0.5), (6, 4, 1, -0.5),
275            // f_234 = 0.5
276            (2, 3, 4, 0.5), (3, 4, 2, 0.5), (4, 2, 3, 0.5),
277            (3, 2, 4, -0.5), (2, 4, 3, -0.5), (4, 3, 2, -0.5),
278            // f_256 = -0.5
279            (2, 5, 6, -0.5), (5, 6, 2, -0.5), (6, 2, 5, -0.5),
280            (5, 2, 6, 0.5), (2, 6, 5, 0.5), (6, 5, 2, 0.5),
281            // f_347 = √3/2
282            (3, 4, 7, SQRT3_HALF), (4, 7, 3, SQRT3_HALF), (7, 3, 4, SQRT3_HALF),
283            (4, 3, 7, -SQRT3_HALF), (3, 7, 4, -SQRT3_HALF), (7, 4, 3, -SQRT3_HALF),
284            // f_567 = √3/2
285            (5, 6, 7, SQRT3_HALF), (6, 7, 5, SQRT3_HALF), (7, 5, 6, SQRT3_HALF),
286            (6, 5, 7, -SQRT3_HALF), (5, 7, 6, -SQRT3_HALF), (7, 6, 5, -SQRT3_HALF),
287        ];
288
289        let mut result = [0.0; 8];
290
291        for &(i, j, k, f) in &STRUCTURE_CONSTANTS {
292            result[k] += self.0[i] * other.0[j] * f;
293        }
294
295        // Apply -2 factor from [X,Y] = -2i Σ f_ijk X_i Y_j T_k
296        for r in &mut result {
297            *r *= -2.0;
298        }
299
300        Self(result)
301    }
302}
303
304// ============================================================================
305// Casimir Operators for SU(3)
306// ============================================================================
307
308impl crate::Casimir for Su3Algebra {
309    type Representation = crate::Su3Irrep;
310
311    fn quadratic_casimir_eigenvalue(irrep: &Self::Representation) -> f64 {
312        let p = irrep.p as f64;
313        let q = irrep.q as f64;
314
315        // c₂(p,q) = (1/3)(p² + q² + pq + 3p + 3q)
316        (p * p + q * q + p * q + 3.0 * p + 3.0 * q) / 3.0
317    }
318
319    /// Cubic Casimir eigenvalue for SU(3).
320    ///
321    /// For representation (p, q), the cubic Casimir eigenvalue is:
322    /// ```text
323    /// c₃(p,q) = (1/18)(p - q)(2p + q + 3)(p + 2q + 3)
324    /// ```
325    ///
326    /// # Properties
327    ///
328    /// - **Conjugation**: c₃(p,q) = -c₃(q,p). Conjugate representations have opposite c₃.
329    /// - **Self-conjugate representations**: c₃ = 0 when p = q (e.g., adjoint (1,1)).
330    /// - **Physical interpretation**: In QCD, distinguishes quarks (positive c₃) from
331    ///   antiquarks (negative c₃) beyond the quadratic Casimir.
332    ///
333    /// # Examples
334    ///
335    /// ```text
336    /// (0,0): c₃ = 0         (trivial)
337    /// (1,0): c₃ = 10/9      (fundamental)
338    /// (0,1): c₃ = -10/9     (antifundamental)
339    /// (1,1): c₃ = 0         (adjoint, self-conjugate)
340    /// (2,0): c₃ = 70/9      (symmetric tensor)
341    /// ```
342    ///
343    /// # Reference
344    ///
345    /// Georgi, "Lie Algebras in Particle Physics" (1999), Chapter 7.
346    fn higher_casimir_eigenvalues(irrep: &Self::Representation) -> Vec<f64> {
347        let p = irrep.p as f64;
348        let q = irrep.q as f64;
349
350        // c₃(p,q) = (1/18)(p - q)(2p + q + 3)(p + 2q + 3)
351        let c3 = (p - q) * (2.0 * p + q + 3.0) * (p + 2.0 * q + 3.0) / 18.0;
352
353        vec![c3]
354    }
355
356    fn rank() -> usize {
357        2 // SU(3) has rank 2 (dimension of Cartan subalgebra)
358    }
359}
360
361impl Su3Algebra {
362    /// Convert algebra element to 3×3 anti-Hermitian matrix
363    ///
364    /// Returns X = i·∑ⱼ aⱼ·λⱼ where λⱼ are Gell-Mann matrices
365    #[must_use]
366    pub fn to_matrix(&self) -> Array2<Complex64> {
367        let [a1, a2, a3, a4, a5, a6, a7, a8] = self.0;
368        let i = Complex64::new(0.0, 1.0);
369        let sqrt3_inv = 1.0 / 3_f64.sqrt();
370
371        // Build i·∑ⱼ aⱼ·λⱼ
372        let mut matrix = Array2::zeros((3, 3));
373
374        // λ₁ = [[0,1,0],[1,0,0],[0,0,0]]
375        matrix[[0, 1]] += i * a1;
376        matrix[[1, 0]] += i * a1;
377
378        // λ₂ = [[0,-i,0],[i,0,0],[0,0,0]]
379        matrix[[0, 1]] += i * Complex64::new(0.0, -a2); // i·(-i·a₂) = a₂
380        matrix[[1, 0]] += i * Complex64::new(0.0, a2); // i·(i·a₂) = -a₂
381
382        // λ₃ = [[1,0,0],[0,-1,0],[0,0,0]]
383        matrix[[0, 0]] += i * a3;
384        matrix[[1, 1]] += -i * a3;
385
386        // λ₄ = [[0,0,1],[0,0,0],[1,0,0]]
387        matrix[[0, 2]] += i * a4;
388        matrix[[2, 0]] += i * a4;
389
390        // λ₅ = [[0,0,-i],[0,0,0],[i,0,0]]
391        matrix[[0, 2]] += i * Complex64::new(0.0, -a5);
392        matrix[[2, 0]] += i * Complex64::new(0.0, a5);
393
394        // λ₆ = [[0,0,0],[0,0,1],[0,1,0]]
395        matrix[[1, 2]] += i * a6;
396        matrix[[2, 1]] += i * a6;
397
398        // λ₇ = [[0,0,0],[0,0,-i],[0,i,0]]
399        matrix[[1, 2]] += i * Complex64::new(0.0, -a7);
400        matrix[[2, 1]] += i * Complex64::new(0.0, a7);
401
402        // λ₈ = (1/√3)·[[1,0,0],[0,1,0],[0,0,-2]]
403        matrix[[0, 0]] += i * a8 * sqrt3_inv;
404        matrix[[1, 1]] += i * a8 * sqrt3_inv;
405        matrix[[2, 2]] += -i * a8 * sqrt3_inv * 2.0;
406
407        matrix
408    }
409
410    /// Extract algebra element from 3×3 anti-Hermitian matrix
411    ///
412    /// Inverse of `to_matrix()`. Uses the normalization Tr(λⱼ λₖ) = 2δⱼₖ.
413    ///
414    /// Given X = i·∑ⱼ aⱼ·λⱼ, we have Tr(X·λⱼ) = i·2·aⱼ, so aⱼ = -i/2·Tr(X·λⱼ).
415    #[must_use]
416    pub fn from_matrix(matrix: &Array2<Complex64>) -> Self {
417        let i = Complex64::new(0.0, 1.0);
418        let neg_i_half = -i / 2.0;
419
420        let mut coeffs = [0.0; 8];
421
422        // For each Gell-Mann matrix, compute aⱼ = -i/2·Tr(X·λⱼ)
423        for j in 0..8 {
424            let lambda_j = Self::gell_mann_matrix(j);
425            let product = matrix.dot(&lambda_j);
426
427            // Compute trace
428            let trace = product[[0, 0]] + product[[1, 1]] + product[[2, 2]];
429
430            // Extract coefficient: aⱼ = -i/2·Tr(X·λⱼ)
431            coeffs[j] = (neg_i_half * trace).re;
432        }
433
434        Self(coeffs)
435    }
436
437    /// Get the j-th Gell-Mann matrix (j = 0..7 for λ₁..λ₈)
438    ///
439    /// Returns the Hermitian Gell-Mann matrix λⱼ (not i·λⱼ).
440    #[must_use]
441    pub fn gell_mann_matrix(j: usize) -> Array2<Complex64> {
442        assert!(j < 8, "Gell-Mann matrices are indexed 0..7");
443
444        let mut matrix = Array2::zeros((3, 3));
445        let i = Complex64::new(0.0, 1.0);
446        let sqrt3_inv = 1.0 / 3_f64.sqrt();
447
448        match j {
449            0 => {
450                // λ₁ = [[0,1,0],[1,0,0],[0,0,0]]
451                matrix[[0, 1]] = Complex64::new(1.0, 0.0);
452                matrix[[1, 0]] = Complex64::new(1.0, 0.0);
453            }
454            1 => {
455                // λ₂ = [[0,-i,0],[i,0,0],[0,0,0]]
456                matrix[[0, 1]] = -i;
457                matrix[[1, 0]] = i;
458            }
459            2 => {
460                // λ₃ = [[1,0,0],[0,-1,0],[0,0,0]]
461                matrix[[0, 0]] = Complex64::new(1.0, 0.0);
462                matrix[[1, 1]] = Complex64::new(-1.0, 0.0);
463            }
464            3 => {
465                // λ₄ = [[0,0,1],[0,0,0],[1,0,0]]
466                matrix[[0, 2]] = Complex64::new(1.0, 0.0);
467                matrix[[2, 0]] = Complex64::new(1.0, 0.0);
468            }
469            4 => {
470                // λ₅ = [[0,0,-i],[0,0,0],[i,0,0]]
471                matrix[[0, 2]] = -i;
472                matrix[[2, 0]] = i;
473            }
474            5 => {
475                // λ₆ = [[0,0,0],[0,0,1],[0,1,0]]
476                matrix[[1, 2]] = Complex64::new(1.0, 0.0);
477                matrix[[2, 1]] = Complex64::new(1.0, 0.0);
478            }
479            6 => {
480                // λ₇ = [[0,0,0],[0,0,-i],[0,i,0]]
481                matrix[[1, 2]] = -i;
482                matrix[[2, 1]] = i;
483            }
484            7 => {
485                // λ₈ = (1/√3)·[[1,0,0],[0,1,0],[0,0,-2]]
486                matrix[[0, 0]] = Complex64::new(sqrt3_inv, 0.0);
487                matrix[[1, 1]] = Complex64::new(sqrt3_inv, 0.0);
488                matrix[[2, 2]] = Complex64::new(-2.0 * sqrt3_inv, 0.0);
489            }
490            _ => unreachable!(),
491        }
492
493        matrix
494    }
495
496    /// SU(3) structure constants `f_ijk`
497    ///
498    /// Returns the structure constant `f_ijk` satisfying [λᵢ, λⱼ] = 2i·∑ₖ fᵢⱼₖ·λₖ
499    /// where λᵢ are the Gell-Mann matrices (0-indexed: 0..7 for λ₁..λ₈).
500    ///
501    /// # Properties
502    ///
503    /// - Totally antisymmetric: `f_ijk` = -`f_jik` = -`f_ikj` = `f_jki`
504    /// - Most entries are zero (only ~24 non-zero out of 512)
505    ///
506    /// # Performance
507    ///
508    /// This is a compile-time constant lookup (O(1)) compared to O(n³) matrix multiplication.
509    ///
510    /// # Note
511    ///
512    /// This function is kept for documentation and testing purposes.
513    /// The `bracket()` method uses a pre-computed table for better performance.
514    #[inline]
515    #[must_use]
516    #[allow(dead_code)]
517    fn structure_constant(i: usize, j: usize, k: usize) -> f64 {
518        const SQRT3_HALF: f64 = 0.866_025_403_784_438_6; // √3/2
519
520        // f_ijk is totally antisymmetric, so f_iik = f_iii = 0
521        if i == j || i == k || j == k {
522            return 0.0;
523        }
524
525        // Use antisymmetry to canonicalize to i < j
526        if i > j {
527            return -Self::structure_constant(j, i, k);
528        }
529
530        // Now i < j, check all non-zero structure constants
531        // (Using 0-based indexing where λ₁→0, λ₂→1, ..., λ₈→7)
532        match (i, j, k) {
533            // f₀₁₂ = 1  ([λ₁, λ₂] = 2i·λ₃)
534            (0, 1, 2) | (1, 2, 0) | (2, 0, 1) => 1.0,
535            (0, 2, 1) | (2, 1, 0) | (1, 0, 2) => -1.0,
536
537            // f₀₃₆ = 1/2  ([λ₁, λ₄] = i·λ₇)
538            (0, 3, 6) | (3, 6, 0) | (6, 0, 3) => 0.5,
539            (0, 6, 3) | (6, 3, 0) | (3, 0, 6) => -0.5,
540
541            // f₀₄₅ = -1/2  ([λ₁, λ₅] = -i·λ₆)
542            (0, 4, 5) | (4, 5, 0) | (5, 0, 4) => -0.5,
543            (0, 5, 4) | (5, 4, 0) | (4, 0, 5) => 0.5,
544
545            // f₁₃₅ = 1/2  ([λ₂, λ₄] = i·λ₆)
546            (1, 3, 5) | (3, 5, 1) | (5, 1, 3) => 0.5,
547            (1, 5, 3) | (5, 3, 1) | (3, 1, 5) => -0.5,
548
549            // f₁₄₆ = 1/2  ([λ₂, λ₅] = i·λ₇)
550            (1, 4, 6) | (4, 6, 1) | (6, 1, 4) => 0.5,
551            (1, 6, 4) | (6, 4, 1) | (4, 1, 6) => -0.5,
552
553            // f₂₃₄ = 1/2  ([λ₃, λ₄] = i·λ₅)
554            (2, 3, 4) | (3, 4, 2) | (4, 2, 3) => 0.5,
555            (2, 4, 3) | (4, 3, 2) | (3, 2, 4) => -0.5,
556
557            // f₂₅₆ = -1/2  ([λ₃, λ₆] = -i·λ₇)
558            (2, 5, 6) | (5, 6, 2) | (6, 2, 5) => -0.5,
559            (2, 6, 5) | (6, 5, 2) | (5, 2, 6) => 0.5,
560
561            // f₃₄₇ = √3/2  ([λ₄, λ₅] = i√3·λ₈)
562            (3, 4, 7) | (4, 7, 3) | (7, 3, 4) => SQRT3_HALF,
563            (3, 7, 4) | (7, 4, 3) | (4, 3, 7) => -SQRT3_HALF,
564
565            // f₅₆₇ = √3/2  ([λ₆, λ₇] = i√3·λ₈)
566            (5, 6, 7) | (6, 7, 5) | (7, 5, 6) => SQRT3_HALF,
567            (5, 7, 6) | (7, 6, 5) | (6, 5, 7) => -SQRT3_HALF,
568
569            // All other combinations are zero
570            _ => 0.0,
571        }
572    }
573}
574
575/// SU(3) group element - 3×3 complex unitary matrix with determinant 1
576///
577/// Represents a color rotation in QCD.
578///
579/// # Representation
580///
581/// We use `Array2<Complex64>` from ndarray to represent the 3×3 unitary matrix.
582///
583/// # Constraints
584///
585/// - Unitarity: U† U = I
586/// - Determinant: det(U) = 1
587///
588/// # Examples
589///
590/// ```
591/// use lie_groups::su3::SU3;
592/// use lie_groups::traits::LieGroup;
593///
594/// let id = SU3::identity();
595/// assert!(id.verify_unitarity(1e-10));
596/// ```
597#[derive(Clone, Debug, PartialEq)]
598pub struct SU3 {
599    /// 3×3 complex unitary matrix
600    pub matrix: Array2<Complex64>,
601}
602
603impl SU3 {
604    /// Identity element
605    #[must_use]
606    pub fn identity() -> Self {
607        let mut matrix = Array2::zeros((3, 3));
608        matrix[[0, 0]] = Complex64::new(1.0, 0.0);
609        matrix[[1, 1]] = Complex64::new(1.0, 0.0);
610        matrix[[2, 2]] = Complex64::new(1.0, 0.0);
611        Self { matrix }
612    }
613
614    /// Verify unitarity: U† U = I
615    #[must_use]
616    pub fn verify_unitarity(&self, tolerance: f64) -> bool {
617        let adjoint = self.matrix.t().mapv(|z| z.conj());
618        let product = adjoint.dot(&self.matrix);
619
620        let mut identity = Array2::zeros((3, 3));
621        identity[[0, 0]] = Complex64::new(1.0, 0.0);
622        identity[[1, 1]] = Complex64::new(1.0, 0.0);
623        identity[[2, 2]] = Complex64::new(1.0, 0.0);
624
625        let diff = product - identity;
626        let norm: f64 = diff
627            .iter()
628            .map(num_complex::Complex::norm_sqr)
629            .sum::<f64>()
630            .sqrt();
631
632        norm < tolerance
633    }
634
635    /// Matrix inverse (equals conjugate transpose for unitary matrices)
636    #[must_use]
637    pub fn inverse(&self) -> Self {
638        Self {
639            matrix: self.matrix.t().mapv(|z| z.conj()),
640        }
641    }
642
643    /// Conjugate transpose (adjoint)
644    #[must_use]
645    pub fn adjoint(&self) -> Self {
646        self.inverse()
647    }
648
649    /// Distance from identity element
650    #[must_use]
651    pub fn distance_to_identity(&self) -> f64 {
652        // Use Frobenius norm: ||U - I||_F
653        let mut identity = Array2::zeros((3, 3));
654        identity[[0, 0]] = Complex64::new(1.0, 0.0);
655        identity[[1, 1]] = Complex64::new(1.0, 0.0);
656        identity[[2, 2]] = Complex64::new(1.0, 0.0);
657
658        let diff = &self.matrix - &identity;
659        diff.iter()
660            .map(num_complex::Complex::norm_sqr)
661            .sum::<f64>()
662            .sqrt()
663    }
664
665    /// Gram-Schmidt reorthogonalization for SU(3) matrices
666    ///
667    /// Projects a potentially corrupted matrix back onto the SU(3) manifold
668    /// using Gram-Schmidt orthogonalization followed by determinant correction.
669    ///
670    /// # Algorithm
671    ///
672    /// 1. Orthogonalize columns using Modified Gram-Schmidt (MGS)
673    /// 2. Normalize to ensure unitarity
674    /// 3. Adjust phase to ensure det(U) = 1
675    ///
676    /// This avoids the log-exp round-trip that would cause infinite recursion
677    /// when called from within `exp()`.
678    ///
679    /// # Numerical Stability
680    ///
681    /// Uses Modified Gram-Schmidt (not Classical GS) for better numerical stability.
682    /// Key difference: projections are computed against already-orthonormalized
683    /// vectors (`result.column(k)`) rather than original columns, and each
684    /// projection is subtracted immediately before the next is computed.
685    /// This provides O(ε) backward error vs O(κε) for Classical GS,
686    /// where κ is the condition number.
687    ///
688    /// # Numerical Considerations (Tao priority)
689    ///
690    /// Uses **scale-relative threshold** for detecting linear dependence:
691    /// ```text
692    /// threshold = max(ε_mach × ||A||_F, ε_abs)
693    /// ```
694    /// where `ε_mach` ≈ 2.2e-16 (machine epsilon) and `ε_abs` = 1e-14 (absolute floor).
695    ///
696    /// This ensures correct behavior for matrices of any scale:
697    /// - For ||A|| ~ 1: threshold ≈ 1e-14 (absolute dominates)
698    /// - For ||A|| ~ 1e-8: threshold ≈ 2e-24 (relative dominates)
699    /// - For ||A|| ~ 1e8: threshold ≈ 2e-8 (relative dominates)
700    ///
701    /// Reference: Björck, "Numerical Methods for Least Squares Problems" (1996)
702    /// Reference: Higham, "Accuracy and Stability of Numerical Algorithms" (2002)
703    #[must_use]
704    fn gram_schmidt_project(matrix: Array2<Complex64>) -> Array2<Complex64> {
705        let mut result: Array2<Complex64> = Array2::zeros((3, 3));
706
707        // Compute Frobenius norm for scale-relative threshold
708        let matrix_norm: f64 = matrix
709            .iter()
710            .map(num_complex::Complex::norm_sqr)
711            .sum::<f64>()
712            .sqrt();
713
714        // Scale-relative threshold: max(ε_mach × ||A||, ε_abs)
715        // This prevents false positives for small-scale matrices and
716        // false negatives for large-scale matrices
717        const MACHINE_EPSILON: f64 = 2.2e-16;
718        const ABSOLUTE_FLOOR: f64 = 1e-14;
719        let relative_threshold = MACHINE_EPSILON * matrix_norm;
720        let threshold = relative_threshold.max(ABSOLUTE_FLOOR);
721
722        // Modified Gram-Schmidt on columns
723        for j in 0..3 {
724            let mut col = matrix.column(j).to_owned();
725
726            // Subtract projections onto previous columns
727            for k in 0..j {
728                let prev_col = result.column(k);
729                let proj: Complex64 = prev_col
730                    .iter()
731                    .zip(col.iter())
732                    .map(|(p, c)| p.conj() * c)
733                    .sum();
734                for i in 0..3 {
735                    col[i] -= proj * prev_col[i];
736                }
737            }
738
739            // Normalize
740            let norm: f64 = col
741                .iter()
742                .map(num_complex::Complex::norm_sqr)
743                .sum::<f64>()
744                .sqrt();
745
746            // Detect linear dependence using scale-relative threshold
747            debug_assert!(
748                norm > threshold,
749                "Gram-Schmidt: column {} is linearly dependent (norm = {:.2e}, threshold = {:.2e}). \
750                 Input matrix is rank-deficient.",
751                j,
752                norm,
753                threshold
754            );
755
756            if norm > threshold {
757                for i in 0..3 {
758                    result[[i, j]] = col[i] / norm;
759                }
760            }
761            // Note: if norm ≤ threshold, column remains zero → det will be ~0 → identity fallback
762        }
763
764        // Ensure det = 1 (SU(N) condition)
765        // Compute determinant and divide by its phase
766        let det = result[[0, 0]]
767            * (result[[1, 1]] * result[[2, 2]] - result[[1, 2]] * result[[2, 1]])
768            - result[[0, 1]] * (result[[1, 0]] * result[[2, 2]] - result[[1, 2]] * result[[2, 0]])
769            + result[[0, 2]] * (result[[1, 0]] * result[[2, 1]] - result[[1, 1]] * result[[2, 0]]);
770
771        // Guard against zero determinant (degenerate matrix)
772        // Use same scale-relative threshold for consistency
773        let det_norm = det.norm();
774        if det_norm < threshold {
775            // Matrix is degenerate; return identity as fallback
776            // This can occur if input was already corrupted
777            return Array2::eye(3);
778        }
779
780        let det_phase = det / det_norm;
781
782        // Multiply by det_phase^{-1/3} to preserve volume while fixing det
783        let correction = (det_phase.conj()).powf(1.0 / 3.0);
784        result.mapv_inplace(|z| z * correction);
785
786        result
787    }
788
789    /// Compute matrix exponential using Taylor series
790    ///
791    /// Assumes ||matrix|| is small (≤ 0.5) for rapid convergence.
792    /// This is a helper method for the scaling-and-squaring algorithm.
793    ///
794    /// # Arguments
795    ///
796    /// - `matrix`: Anti-Hermitian 3×3 matrix
797    /// - `max_terms`: Maximum number of Taylor series terms
798    ///
799    /// # Returns
800    ///
801    /// exp(matrix) ∈ SU(3)
802    #[must_use]
803    fn exp_taylor(matrix: &Array2<Complex64>, max_terms: usize) -> Self {
804        let mut result = Array2::zeros((3, 3));
805        result[[0, 0]] = Complex64::new(1.0, 0.0);
806        result[[1, 1]] = Complex64::new(1.0, 0.0);
807        result[[2, 2]] = Complex64::new(1.0, 0.0);
808
809        let mut term = matrix.clone();
810        let mut factorial = 1.0;
811
812        for n in 1..=max_terms {
813            factorial *= n as f64;
814            result += &term.mapv(|z| z / factorial);
815
816            // Early termination if term is negligible
817            let term_norm: f64 = term
818                .iter()
819                .map(num_complex::Complex::norm_sqr)
820                .sum::<f64>()
821                .sqrt();
822
823            if term_norm / factorial < 1e-14 {
824                break;
825            }
826
827            if n < max_terms {
828                term = term.dot(matrix);
829            }
830        }
831
832        Self { matrix: result }
833    }
834
835    /// Compute matrix square root using Denman-Beavers iteration.
836    ///
837    /// For a unitary matrix U, this converges to U^{1/2} (principal square root).
838    ///
839    /// # Algorithm
840    ///
841    /// Y₀ = U, Z₀ = I
842    /// Yₙ₊₁ = (Yₙ + Zₙ⁻¹)/2
843    /// Zₙ₊₁ = (Zₙ + Yₙ⁻¹)/2
844    ///
845    /// Converges quadratically to Y = U^{1/2}, Z = U^{-1/2}.
846    ///
847    /// # Reference
848    ///
849    /// Higham, "Functions of Matrices", Ch. 6.
850    fn matrix_sqrt_db(u: &Array2<Complex64>) -> Array2<Complex64> {
851        use nalgebra::{Complex as NaComplex, Matrix3};
852
853        // Convert ndarray to nalgebra for inversion
854        fn to_nalgebra(a: &Array2<Complex64>) -> Matrix3<NaComplex<f64>> {
855            Matrix3::from_fn(|i, j| NaComplex::new(a[[i, j]].re, a[[i, j]].im))
856        }
857
858        fn to_ndarray(m: &Matrix3<NaComplex<f64>>) -> Array2<Complex64> {
859            Array2::from_shape_fn((3, 3), |(i, j)| Complex64::new(m[(i, j)].re, m[(i, j)].im))
860        }
861
862        let mut y = to_nalgebra(u);
863        let mut z = Matrix3::<NaComplex<f64>>::identity();
864
865        const MAX_ITERS: usize = 20;
866        const TOL: f64 = 1e-14;
867
868        for _ in 0..MAX_ITERS {
869            // Compute actual inverses
870            let y_inv = y.try_inverse().unwrap_or(y.adjoint()); // Fallback to adjoint if singular
871            let z_inv = z.try_inverse().unwrap_or(z.adjoint());
872
873            let y_new = (y + z_inv).scale(0.5);
874            let z_new = (z + y_inv).scale(0.5);
875
876            // Check convergence
877            let diff: f64 = (y_new - y).norm();
878
879            y = y_new;
880            z = z_new;
881
882            if diff < TOL {
883                break;
884            }
885        }
886
887        to_ndarray(&y)
888    }
889}
890
891impl fmt::Display for Su3Algebra {
892    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
893        write!(f, "su(3)[")?;
894        for (i, c) in self.0.iter().enumerate() {
895            if i > 0 {
896                write!(f, ", ")?;
897            }
898            write!(f, "{:.4}", c)?;
899        }
900        write!(f, "]")
901    }
902}
903
904impl fmt::Display for SU3 {
905    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
906        let dist = self.distance_to_identity();
907        write!(f, "SU(3)(d={:.4})", dist)
908    }
909}
910
911/// Group multiplication: U₁ · U₂
912impl Mul<&SU3> for &SU3 {
913    type Output = SU3;
914    fn mul(self, rhs: &SU3) -> SU3 {
915        SU3 {
916            matrix: self.matrix.dot(&rhs.matrix),
917        }
918    }
919}
920
921impl Mul<&SU3> for SU3 {
922    type Output = SU3;
923    fn mul(self, rhs: &SU3) -> SU3 {
924        &self * rhs
925    }
926}
927
928impl MulAssign<&SU3> for SU3 {
929    fn mul_assign(&mut self, rhs: &SU3) {
930        self.matrix = self.matrix.dot(&rhs.matrix);
931    }
932}
933
934impl LieGroup for SU3 {
935    const DIM: usize = 3;
936
937    type Algebra = Su3Algebra;
938
939    fn identity() -> Self {
940        Self::identity()
941    }
942
943    fn compose(&self, other: &Self) -> Self {
944        Self {
945            matrix: self.matrix.dot(&other.matrix),
946        }
947    }
948
949    fn inverse(&self) -> Self {
950        Self::inverse(self)
951    }
952
953    fn adjoint(&self) -> Self {
954        Self::adjoint(self)
955    }
956
957    fn adjoint_action(&self, algebra_element: &Su3Algebra) -> Su3Algebra {
958        // Ad_g(X) = g X g† for matrix groups
959        let x_matrix = algebra_element.to_matrix();
960        let g_x = self.matrix.dot(&x_matrix);
961        let g_adjoint_matrix = self.matrix.t().mapv(|z| z.conj());
962        let result = g_x.dot(&g_adjoint_matrix);
963
964        Su3Algebra::from_matrix(&result)
965    }
966
967    fn distance_to_identity(&self) -> f64 {
968        Self::distance_to_identity(self)
969    }
970
971    fn exp(tangent: &Su3Algebra) -> Self {
972        // Matrix exponential using scaling-and-squaring algorithm
973        //
974        // Algorithm: exp(X) = [exp(X/2^k)]^(2^k)
975        // where k is chosen so ||X/2^k|| is small enough for Taylor series convergence
976        //
977        // This ensures accurate computation even for large ||X||
978
979        let x_matrix = tangent.to_matrix();
980
981        // Compute matrix norm (Frobenius norm)
982        let norm: f64 = x_matrix
983            .iter()
984            .map(num_complex::Complex::norm_sqr)
985            .sum::<f64>()
986            .sqrt();
987
988        // Choose scaling parameter k such that ||X/2^k|| ≤ 0.5
989        // This ensures rapid Taylor series convergence
990        let k = if norm > 0.5 {
991            (norm / 0.5).log2().ceil() as usize
992        } else {
993            0
994        };
995
996        // Scale the matrix: Y = X / 2^k
997        let scale_factor = 1.0 / (1_u64 << k) as f64;
998        let scaled_matrix = x_matrix.mapv(|z| z * scale_factor);
999
1000        // Compute exp(Y) using Taylor series (converges rapidly for ||Y|| ≤ 0.5)
1001        let exp_scaled = SU3::exp_taylor(&scaled_matrix, 15);
1002
1003        // Square k times to recover exp(X) = [exp(Y)]^(2^k)
1004        //
1005        // Numerical stability (Tao priority): Reorthogonalize after EVERY squaring.
1006        //
1007        // Rationale (Higham & Al-Mohy, 2010):
1008        // - Each matrix multiplication accumulates O(nε) orthogonality loss
1009        // - After k squarings without reorthogonalization: O(2^k × nε) error
1010        // - For k=10, this is O(1000ε) ≈ 1e-13, approaching catastrophic loss
1011        //
1012        // Previous code reorthogonalized every 4 squarings, which is insufficient
1013        // for large k or ill-conditioned intermediate results. The cost of
1014        // Gram-Schmidt (O(n³) = O(27) for 3×3) is negligible compared to
1015        // the matrix multiply (also O(n³)).
1016        let mut result = exp_scaled.matrix;
1017        for _ in 0..k {
1018            result = result.dot(&result);
1019            // Reorthogonalize after every squaring to maintain manifold constraint
1020            result = Self::gram_schmidt_project(result);
1021        }
1022
1023        // Result is already orthogonalized from final loop iteration
1024        Self { matrix: result }
1025    }
1026
1027    fn log(&self) -> crate::error::LogResult<Su3Algebra> {
1028        use crate::error::LogError;
1029
1030        // Matrix logarithm for SU(3) using inverse scaling-squaring algorithm.
1031        //
1032        // Algorithm (Higham, "Functions of Matrices", Ch. 11):
1033        // 1. Take square roots until ||U^{1/2^k} - I|| < 0.5
1034        // 2. Use Taylor series for log(I + X) with ||X|| < 0.5 (fast convergence)
1035        // 3. Scale back: log(U) = 2^k × log(U^{1/2^k})
1036        //
1037        // This achieves ~1e-10 accuracy vs ~1e-2 for direct Taylor series.
1038
1039        // Check distance from identity
1040        let dist = self.distance_to_identity();
1041        const MAX_DISTANCE: f64 = 2.0; // Maximum distance for principal branch
1042
1043        if dist > MAX_DISTANCE {
1044            return Err(LogError::NotNearIdentity {
1045                distance: dist,
1046                threshold: MAX_DISTANCE,
1047            });
1048        }
1049
1050        // Check if at identity
1051        if dist < 1e-14 {
1052            return Ok(Su3Algebra::zero());
1053        }
1054
1055        // Phase 1: Inverse scaling via matrix square roots
1056        // Take square roots until ||U - I|| < 0.5 for rapid Taylor convergence
1057        let mut current = self.matrix.clone();
1058        let mut num_sqrts = 0;
1059        const MAX_SQRTS: usize = 32; // Prevent infinite loop
1060        const TARGET_NORM: f64 = 0.5; // Taylor converges rapidly for ||X|| < 0.5
1061
1062        let identity: Array2<Complex64> = Array2::eye(3);
1063
1064        while num_sqrts < MAX_SQRTS {
1065            let x_matrix = &current - &identity;
1066            let x_norm: f64 = x_matrix
1067                .iter()
1068                .map(num_complex::Complex::norm_sqr)
1069                .sum::<f64>()
1070                .sqrt();
1071
1072            if x_norm < TARGET_NORM {
1073                break;
1074            }
1075
1076            // Compute matrix square root using Denman-Beavers iteration
1077            current = Self::matrix_sqrt_db(&current);
1078            num_sqrts += 1;
1079        }
1080
1081        // Phase 2: Taylor series for log(I + X) with ||X|| < 0.5
1082        let x_matrix = &current - &identity;
1083
1084        // Taylor series: log(I + X) = Σ_{n=1}^∞ (-1)^{n+1} X^n / n
1085        let mut log_matrix = x_matrix.clone();
1086        let mut x_power = x_matrix.clone();
1087
1088        // With ||X|| < 0.5, 30 terms gives ~0.5^30/30 ≈ 3e-11 truncation error
1089        const N_TERMS: usize = 30;
1090
1091        for n in 2..=N_TERMS {
1092            x_power = x_power.dot(&x_matrix);
1093            let coefficient = (-1.0_f64).powi(n as i32 + 1) / n as f64;
1094            log_matrix = log_matrix + x_power.mapv(|z| z * coefficient);
1095        }
1096
1097        // Phase 3: Scale back: log(U) = 2^k × log(U^{1/2^k})
1098        let scale_factor = (1_u64 << num_sqrts) as f64;
1099        log_matrix = log_matrix.mapv(|z| z * scale_factor);
1100
1101        // Convert result to algebra element
1102        Ok(Su3Algebra::from_matrix(&log_matrix))
1103    }
1104
1105    fn dim() -> usize {
1106        3 // SU(3) consists of 3×3 matrices
1107    }
1108
1109    fn trace(&self) -> Complex64 {
1110        self.matrix[[0, 0]] + self.matrix[[1, 1]] + self.matrix[[2, 2]]
1111    }
1112}
1113
1114// ============================================================================
1115// Mathematical Property Implementations
1116// ============================================================================
1117
1118use crate::traits::{Compact, SemiSimple, Simple};
1119
1120/// SU(3) is compact
1121///
1122/// All elements are bounded: ||U|| = 1 for all U ∈ SU(3).
1123impl Compact for SU3 {}
1124
1125/// SU(3) is simple
1126///
1127/// It has no non-trivial normal subgroups (except center ℤ₃).
1128impl Simple for SU3 {}
1129
1130/// SU(3) is semi-simple
1131impl SemiSimple for SU3 {}
1132
1133// ============================================================================
1134// Algebra Marker Traits
1135// ============================================================================
1136
1137/// su(3) algebra elements are traceless by construction.
1138///
1139/// The representation `Su3Algebra([f64; 8])` stores coefficients in the
1140/// Gell-Mann basis {iλ₁, ..., iλ₈}. All Gell-Mann matrices are traceless.
1141impl TracelessByConstruction for Su3Algebra {}
1142
1143/// su(3) algebra elements are anti-Hermitian by construction.
1144///
1145/// The representation uses {iλⱼ} where λⱼ are Hermitian Gell-Mann matrices.
1146impl AntiHermitianByConstruction for Su3Algebra {}
1147
1148#[cfg(test)]
1149mod tests {
1150    use super::*;
1151    use approx::assert_relative_eq;
1152
1153    #[test]
1154    fn test_identity() {
1155        let id = SU3::identity();
1156        assert!(id.verify_unitarity(1e-10));
1157        assert_relative_eq!(id.distance_to_identity(), 0.0, epsilon = 1e-10);
1158    }
1159
1160    #[test]
1161    fn test_algebra_dimension() {
1162        assert_eq!(Su3Algebra::dim(), 8);
1163    }
1164
1165    #[test]
1166    fn test_gell_mann_hermiticity() {
1167        // Verify that Gell-Mann matrices are Hermitian: λ† = λ
1168        for j in 0..8 {
1169            let lambda = Su3Algebra::gell_mann_matrix(j);
1170            let adjoint = lambda.t().mapv(|z| z.conj());
1171
1172            for i in 0..3 {
1173                for k in 0..3 {
1174                    assert_relative_eq!(lambda[[i, k]].re, adjoint[[i, k]].re, epsilon = 1e-10);
1175                    assert_relative_eq!(lambda[[i, k]].im, adjoint[[i, k]].im, epsilon = 1e-10);
1176                }
1177            }
1178        }
1179    }
1180
1181    #[test]
1182    fn test_gell_mann_traceless() {
1183        // Verify that Gell-Mann matrices are traceless
1184        for j in 0..8 {
1185            let lambda = Su3Algebra::gell_mann_matrix(j);
1186            let trace = lambda[[0, 0]] + lambda[[1, 1]] + lambda[[2, 2]];
1187            assert_relative_eq!(trace.re, 0.0, epsilon = 1e-10);
1188            assert_relative_eq!(trace.im, 0.0, epsilon = 1e-10);
1189        }
1190    }
1191
1192    #[test]
1193    fn test_matrix_roundtrip() {
1194        // Test that to_matrix() and from_matrix() are inverses
1195        let algebra = Su3Algebra([1.0, 2.0, 3.0, 0.5, -0.5, 1.5, -1.5, 0.3]);
1196        let matrix = algebra.to_matrix();
1197        let recovered = Su3Algebra::from_matrix(&matrix);
1198
1199        for i in 0..8 {
1200            assert_relative_eq!(algebra.0[i], recovered.0[i], epsilon = 1e-10);
1201        }
1202    }
1203
1204    #[test]
1205    fn test_inverse() {
1206        use crate::traits::LieGroup;
1207
1208        let g = SU3::exp(&Su3Algebra([0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8]));
1209        let g_inv = g.inverse();
1210        let product = g.compose(&g_inv);
1211
1212        assert_relative_eq!(product.distance_to_identity(), 0.0, epsilon = 1e-6);
1213    }
1214
1215    #[test]
1216    fn test_adjoint_identity() {
1217        use crate::traits::LieGroup;
1218
1219        let e = SU3::identity();
1220        let x = Su3Algebra([1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]);
1221        let result = e.adjoint_action(&x);
1222
1223        for i in 0..8 {
1224            assert_relative_eq!(result.0[i], x.0[i], epsilon = 1e-10);
1225        }
1226    }
1227
1228    #[test]
1229    fn test_structure_constants_bracket() {
1230        use crate::traits::LieAlgebra;
1231
1232        // Test [λ₁, λ₂] = 2i·λ₃ ⟹ coefficients: [λ₁, λ₂] has c₃ ≠ 0
1233        let lambda1 = Su3Algebra::basis_element(0);
1234        let lambda2 = Su3Algebra::basis_element(1);
1235        let bracket = lambda1.bracket(&lambda2);
1236
1237        // Should get result with λ₃ component (index 2)
1238        assert_relative_eq!(bracket.0[2], -2.0, epsilon = 1e-10); // Factor of -2 from formula
1239        for i in [0, 1, 3, 4, 5, 6, 7] {
1240            assert_relative_eq!(bracket.0[i], 0.0, epsilon = 1e-10);
1241        }
1242
1243        // Test antisymmetry: [λ₂, λ₁] = -[λ₁, λ₂]
1244        let bracket_reversed = lambda2.bracket(&lambda1);
1245        for i in 0..8 {
1246            assert_relative_eq!(bracket.0[i], -bracket_reversed.0[i], epsilon = 1e-10);
1247        }
1248
1249        // Test [λ₄, λ₅] = 2i·(√3/2)·λ₈ ⟹ coefficient c₈ = -√3
1250        let lambda4 = Su3Algebra::basis_element(3);
1251        let lambda5 = Su3Algebra::basis_element(4);
1252        let bracket_45 = lambda4.bracket(&lambda5);
1253        let expected_c8 = -2.0 * (3.0_f64.sqrt() / 2.0);
1254        assert_relative_eq!(bracket_45.0[7], expected_c8, epsilon = 1e-10);
1255    }
1256
1257    #[test]
1258    fn test_bracket_jacobi_identity() {
1259        use crate::traits::LieAlgebra;
1260
1261        let x = Su3Algebra::basis_element(0);
1262        let y = Su3Algebra::basis_element(3);
1263        let z = Su3Algebra::basis_element(7);
1264
1265        // [X, [Y, Z]] + [Y, [Z, X]] + [Z, [X, Y]] = 0
1266        let t1 = x.bracket(&y.bracket(&z));
1267        let t2 = y.bracket(&z.bracket(&x));
1268        let t3 = z.bracket(&x.bracket(&y));
1269        let sum = t1.add(&t2).add(&t3);
1270
1271        assert!(
1272            sum.norm() < 1e-10,
1273            "Jacobi identity violated for SU(3): ||sum|| = {:.2e}",
1274            sum.norm()
1275        );
1276    }
1277
1278    #[test]
1279    fn test_bracket_bilinearity() {
1280        use crate::traits::LieAlgebra;
1281
1282        let x = Su3Algebra::basis_element(0);
1283        let y = Su3Algebra::basis_element(2);
1284        let z = Su3Algebra::basis_element(5);
1285        let alpha = 3.7;
1286
1287        // [αX + Y, Z] = α[X, Z] + [Y, Z]
1288        let lhs = x.scale(alpha).add(&y).bracket(&z);
1289        let rhs = x.bracket(&z).scale(alpha).add(&y.bracket(&z));
1290        for i in 0..8 {
1291            assert_relative_eq!(lhs.0[i], rhs.0[i], epsilon = 1e-12);
1292        }
1293    }
1294
1295    #[test]
1296    fn test_exp_large_algebra_element() {
1297        use crate::traits::LieGroup;
1298
1299        // Test with large algebra element (||X|| > 1)
1300        // Old Taylor series would fail, scaling-and-squaring handles this
1301        let large_algebra = Su3Algebra([2.0, 1.5, -1.8, 0.9, -1.2, 1.1, -0.8, 1.3]);
1302        let norm = large_algebra.norm();
1303
1304        // Verify this is actually large
1305        assert!(norm > 1.0, "Test requires ||X|| > 1, got {}", norm);
1306
1307        // Compute exponential (should not panic or produce NaN)
1308        let g = SU3::exp(&large_algebra);
1309
1310        // Verify unitarity is preserved
1311        assert!(
1312            g.verify_unitarity(1e-8),
1313            "Unitarity violated for large algebra element"
1314        );
1315
1316        // Verify g is not identity (non-trivial rotation)
1317        assert!(g.distance_to_identity() > 0.1);
1318    }
1319
1320    #[test]
1321    fn test_exp_very_small_algebra_element() {
1322        use crate::traits::LieGroup;
1323
1324        // Test with very small algebra element
1325        let small_algebra = Su3Algebra([1e-8, 2e-8, -1e-8, 3e-9, -2e-9, 1e-9, -5e-10, 2e-10]);
1326
1327        let g = SU3::exp(&small_algebra);
1328
1329        // Should be very close to identity
1330        assert!(g.distance_to_identity() < 1e-7);
1331        assert!(g.verify_unitarity(1e-12));
1332    }
1333
1334    #[test]
1335    fn test_exp_scaling_correctness() {
1336        use crate::traits::LieGroup;
1337
1338        // Verify exp(2X) = exp(X)^2 (approximately)
1339        let algebra = Su3Algebra([0.5, 0.3, -0.4, 0.2, -0.3, 0.1, -0.2, 0.25]);
1340
1341        let exp_x = SU3::exp(&algebra);
1342        let exp_2x = SU3::exp(&algebra.scale(2.0));
1343        let exp_x_squared = exp_x.compose(&exp_x);
1344
1345        // These should be approximately equal
1346        let distance = exp_2x.distance(&exp_x_squared);
1347        assert!(
1348            distance < 1e-6,
1349            "exp(2X) should equal exp(X)^2, distance = {}",
1350            distance
1351        );
1352    }
1353
1354    // ========================================================================
1355    // Property-Based Tests for Group Axioms
1356    // ========================================================================
1357    //
1358    // These tests use proptest to verify that SU(3) satisfies the
1359    // mathematical axioms of a Lie group for randomly generated elements.
1360    //
1361    // Run with: cargo test --features nightly
1362
1363    #[cfg(feature = "proptest")]
1364    use proptest::prelude::*;
1365
1366    /// Strategy for generating arbitrary SU(3) elements.
1367    ///
1368    /// We generate SU(3) elements via the exponential map from random
1369    /// algebra elements. Using smaller range (-0.5..0.5) for better
1370    /// convergence of the Taylor series exponential.
1371    #[cfg(feature = "proptest")]
1372    fn arb_su3() -> impl Strategy<Value = SU3> {
1373        // Use smaller range for Taylor series convergence
1374        let range = -0.5_f64..0.5_f64;
1375
1376        (
1377            range.clone(),
1378            range.clone(),
1379            range.clone(),
1380            range.clone(),
1381            range.clone(),
1382            range.clone(),
1383            range.clone(),
1384            range,
1385        )
1386            .prop_map(|(a1, a2, a3, a4, a5, a6, a7, a8)| {
1387                let algebra = Su3Algebra([a1, a2, a3, a4, a5, a6, a7, a8]);
1388                SU3::exp(&algebra)
1389            })
1390    }
1391
1392    /// Strategy for generating arbitrary `Su3Algebra` elements.
1393    #[cfg(feature = "proptest")]
1394    fn arb_su3_algebra() -> impl Strategy<Value = Su3Algebra> {
1395        let range = -0.5_f64..0.5_f64;
1396
1397        (
1398            range.clone(),
1399            range.clone(),
1400            range.clone(),
1401            range.clone(),
1402            range.clone(),
1403            range.clone(),
1404            range.clone(),
1405            range,
1406        )
1407            .prop_map(|(a1, a2, a3, a4, a5, a6, a7, a8)| {
1408                Su3Algebra([a1, a2, a3, a4, a5, a6, a7, a8])
1409            })
1410    }
1411
1412    #[cfg(feature = "proptest")]
1413    proptest! {
1414        /// **Group Axiom 1: Identity Element**
1415        ///
1416        /// For all U ∈ SU(3):
1417        /// - I · U = U (left identity)
1418        /// - U · I = U (right identity)
1419        ///
1420        /// where I = identity element
1421        #[test]
1422        fn prop_identity_axiom(u in arb_su3()) {
1423            let e = SU3::identity();
1424
1425            // Left identity: I · U = U
1426            let left = e.compose(&u);
1427            prop_assert!(
1428                left.distance(&u) < 1e-6,
1429                "Left identity failed: I·U != U, distance = {}",
1430                left.distance(&u)
1431            );
1432
1433            // Right identity: U · I = U
1434            let right = u.compose(&e);
1435            prop_assert!(
1436                right.distance(&u) < 1e-6,
1437                "Right identity failed: U·I != U, distance = {}",
1438                right.distance(&u)
1439            );
1440        }
1441
1442        /// **Group Axiom 2: Inverse Element**
1443        ///
1444        /// For all U ∈ SU(3):
1445        /// - U · U† = I (right inverse)
1446        /// - U† · U = I (left inverse)
1447        ///
1448        /// where U† = conjugate transpose
1449        #[test]
1450        fn prop_inverse_axiom(u in arb_su3()) {
1451            let u_inv = u.inverse();
1452
1453            // Right inverse: U · U† = I
1454            let right_product = u.compose(&u_inv);
1455            prop_assert!(
1456                right_product.is_near_identity(1e-6),
1457                "Right inverse failed: U·U† != I, distance = {}",
1458                right_product.distance_to_identity()
1459            );
1460
1461            // Left inverse: U† · U = I
1462            let left_product = u_inv.compose(&u);
1463            prop_assert!(
1464                left_product.is_near_identity(1e-6),
1465                "Left inverse failed: U†·U != I, distance = {}",
1466                left_product.distance_to_identity()
1467            );
1468        }
1469
1470        /// **Group Axiom 3: Associativity**
1471        ///
1472        /// For all U₁, U₂, U₃ ∈ SU(3):
1473        /// - (U₁ · U₂) · U₃ = U₁ · (U₂ · U₃)
1474        ///
1475        /// Group multiplication is associative.
1476        #[test]
1477        fn prop_associativity(u1 in arb_su3(), u2 in arb_su3(), u3 in arb_su3()) {
1478            // Left association: (U₁ · U₂) · U₃
1479            let left_assoc = u1.compose(&u2).compose(&u3);
1480
1481            // Right association: U₁ · (U₂ · U₃)
1482            let right_assoc = u1.compose(&u2.compose(&u3));
1483
1484            prop_assert!(
1485                left_assoc.distance(&right_assoc) < 1e-6,
1486                "Associativity failed: (U₁·U₂)·U₃ != U₁·(U₂·U₃), distance = {}",
1487                left_assoc.distance(&right_assoc)
1488            );
1489        }
1490
1491        /// **Lie Group Property: Inverse is Smooth**
1492        ///
1493        /// For SU(3), the inverse operation is smooth (continuously differentiable).
1494        /// We verify this by checking that nearby elements have nearby inverses.
1495        #[test]
1496        fn prop_inverse_continuity(u in arb_su3()) {
1497            // Create a small perturbation
1498            let epsilon = 0.01;
1499            let perturbation = SU3::exp(&Su3Algebra([epsilon, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]));
1500            let u_perturbed = u.compose(&perturbation);
1501
1502            // Check that inverses are close
1503            let inv_distance = u.inverse().distance(&u_perturbed.inverse());
1504
1505            prop_assert!(
1506                inv_distance < 0.1,
1507                "Inverse not continuous: small perturbation caused large inverse change, distance = {}",
1508                inv_distance
1509            );
1510        }
1511
1512        /// **Unitarity Preservation**
1513        ///
1514        /// All SU(3) operations should preserve unitarity.
1515        /// This is not strictly a group axiom, but it's essential for SU(3).
1516        #[test]
1517        fn prop_unitarity_preserved(u1 in arb_su3(), u2 in arb_su3()) {
1518            // Composition preserves unitarity
1519            let product = u1.compose(&u2);
1520            prop_assert!(
1521                product.verify_unitarity(1e-10),
1522                "Composition violated unitarity"
1523            );
1524
1525            // Inverse preserves unitarity
1526            let inv = u1.inverse();
1527            prop_assert!(
1528                inv.verify_unitarity(1e-10),
1529                "Inverse violated unitarity"
1530            );
1531        }
1532
1533        /// **Adjoint Representation: Group Homomorphism**
1534        ///
1535        /// The adjoint representation Ad: G → Aut(𝔤) is a group homomorphism:
1536        /// - Ad_{U₁∘U₂}(X) = Ad_{U₁}(Ad_{U₂}(X))
1537        ///
1538        /// This is a fundamental property that must hold for the adjoint action
1539        /// to be a valid representation of the group.
1540        #[test]
1541        fn prop_adjoint_homomorphism(
1542            u1 in arb_su3(),
1543            u2 in arb_su3(),
1544            x in arb_su3_algebra()
1545        ) {
1546            // Compute Ad_{U₁∘U₂}(X)
1547            let u_composed = u1.compose(&u2);
1548            let left = u_composed.adjoint_action(&x);
1549
1550            // Compute Ad_{U₁}(Ad_{U₂}(X))
1551            let ad_u2_x = u2.adjoint_action(&x);
1552            let right = u1.adjoint_action(&ad_u2_x);
1553
1554            // They should be equal
1555            let diff = left.add(&right.scale(-1.0));
1556            prop_assert!(
1557                diff.norm() < 1e-6,
1558                "Adjoint homomorphism failed: Ad_{{U₁∘U₂}}(X) != Ad_{{U₁}}(Ad_{{U₂}}(X)), diff norm = {}",
1559                diff.norm()
1560            );
1561        }
1562
1563        /// **Adjoint Representation: Identity Action**
1564        ///
1565        /// The identity element acts trivially on the Lie algebra:
1566        /// - Ad_I(X) = X for all X ∈ su(3)
1567        #[test]
1568        fn prop_adjoint_identity(x in arb_su3_algebra()) {
1569            let e = SU3::identity();
1570            let result = e.adjoint_action(&x);
1571
1572            let diff = result.add(&x.scale(-1.0));
1573            prop_assert!(
1574                diff.norm() < 1e-10,
1575                "Identity action failed: Ad_I(X) != X, diff norm = {}",
1576                diff.norm()
1577            );
1578        }
1579
1580        /// **Adjoint Representation: Lie Bracket Preservation**
1581        ///
1582        /// The adjoint representation preserves the Lie bracket:
1583        /// - Ad_U([X,Y]) = [Ad_U(X), Ad_U(Y)]
1584        ///
1585        /// This is a critical property that ensures the adjoint action
1586        /// is a Lie algebra automorphism.
1587        #[test]
1588        fn prop_adjoint_bracket_preservation(
1589            u in arb_su3(),
1590            x in arb_su3_algebra(),
1591            y in arb_su3_algebra()
1592        ) {
1593            use crate::traits::LieAlgebra;
1594
1595            // Compute Ad_U([X,Y])
1596            let bracket_xy = x.bracket(&y);
1597            let left = u.adjoint_action(&bracket_xy);
1598
1599            // Compute [Ad_U(X), Ad_U(Y)]
1600            let ad_x = u.adjoint_action(&x);
1601            let ad_y = u.adjoint_action(&y);
1602            let right = ad_x.bracket(&ad_y);
1603
1604            // They should be equal
1605            let diff = left.add(&right.scale(-1.0));
1606            prop_assert!(
1607                diff.norm() < 1e-5,
1608                "Bracket preservation failed: Ad_U([X,Y]) != [Ad_U(X), Ad_U(Y)], diff norm = {}",
1609                diff.norm()
1610            );
1611        }
1612
1613        /// **Adjoint Representation: Inverse Property**
1614        ///
1615        /// The inverse of an element acts as the inverse transformation:
1616        /// - Ad_{U†}(Ad_U(X)) = X
1617        #[test]
1618        fn prop_adjoint_inverse(u in arb_su3(), x in arb_su3_algebra()) {
1619            // Apply Ad_U then Ad_{U†}
1620            let ad_u_x = u.adjoint_action(&x);
1621            let u_inv = u.inverse();
1622            let result = u_inv.adjoint_action(&ad_u_x);
1623
1624            // Should recover X
1625            let diff = result.add(&x.scale(-1.0));
1626            prop_assert!(
1627                diff.norm() < 1e-6,
1628                "Inverse property failed: Ad_{{U†}}(Ad_U(X)) != X, diff norm = {}",
1629                diff.norm()
1630            );
1631        }
1632    }
1633
1634    // ========================================================================
1635    // Exp/Log Round-Trip Tests
1636    // ========================================================================
1637
1638    /// **Exp-Log Round-Trip Test for SU(3)**
1639    ///
1640    /// For any algebra element X with small norm, we should have:
1641    /// ```text
1642    /// log(exp(X)) ≈ X
1643    /// ```
1644    #[test]
1645    fn test_exp_log_roundtrip() {
1646        use crate::traits::{LieAlgebra, LieGroup};
1647        use rand::SeedableRng;
1648        use rand_distr::{Distribution, Uniform};
1649
1650        let mut rng = rand::rngs::StdRng::seed_from_u64(11111);
1651        let dist = Uniform::new(-0.5, 0.5); // Small values for stable exp/log
1652
1653        for _ in 0..50 {
1654            let mut coeffs = [0.0; 8];
1655            for coeff in &mut coeffs {
1656                *coeff = dist.sample(&mut rng);
1657            }
1658            let x = Su3Algebra(coeffs);
1659
1660            // exp then log
1661            let g = SU3::exp(&x);
1662            let x_recovered = g.log().expect("log should succeed for exp output");
1663
1664            // Inverse scaling-squaring algorithm achieves machine precision (~1e-14)
1665            let diff = x.add(&x_recovered.scale(-1.0));
1666            assert!(
1667                diff.norm() < 1e-10,
1668                "log(exp(X)) should equal X: ||diff|| = {:.2e}",
1669                diff.norm()
1670            );
1671        }
1672    }
1673
1674    /// **Log-Exp Round-Trip Test for SU(3)**
1675    ///
1676    /// For any group element g, we should have:
1677    /// ```text
1678    /// exp(log(g)) = g
1679    /// ```
1680    #[test]
1681    fn test_log_exp_roundtrip() {
1682        use crate::traits::LieGroup;
1683        use rand::SeedableRng;
1684        use rand_distr::{Distribution, Uniform};
1685
1686        let mut rng = rand::rngs::StdRng::seed_from_u64(22222);
1687        let dist = Uniform::new(-0.5, 0.5);
1688
1689        for _ in 0..50 {
1690            // Generate random SU(3) element via exp of small algebra element
1691            let mut coeffs = [0.0; 8];
1692            for coeff in &mut coeffs {
1693                *coeff = dist.sample(&mut rng);
1694            }
1695            let g = SU3::exp(&Su3Algebra(coeffs));
1696
1697            // log then exp
1698            let x = g.log().expect("log should succeed for valid SU(3) element");
1699            let g_recovered = SU3::exp(&x);
1700
1701            // Inverse scaling-squaring algorithm achieves machine precision (~1e-14)
1702            let diff = g.compose(&g_recovered.inverse()).distance_to_identity();
1703            assert!(
1704                diff < 1e-10,
1705                "exp(log(g)) should equal g: diff = {:.2e}",
1706                diff
1707            );
1708        }
1709    }
1710
1711    // ========================================================================
1712    // Casimir Operator Tests
1713    // ========================================================================
1714
1715    #[test]
1716    fn test_su3_casimir_trivial() {
1717        // (0,0) trivial: c₂ = 0
1718        use crate::Casimir;
1719        use crate::Su3Irrep;
1720
1721        let c2 = Su3Algebra::quadratic_casimir_eigenvalue(&Su3Irrep::TRIVIAL);
1722        assert_eq!(c2, 0.0, "Casimir of trivial representation should be 0");
1723    }
1724
1725    #[test]
1726    fn test_su3_casimir_fundamental() {
1727        // (1,0) fundamental (quark): c₂ = 4/3
1728        use crate::Casimir;
1729        use crate::Su3Irrep;
1730
1731        let c2 = Su3Algebra::quadratic_casimir_eigenvalue(&Su3Irrep::FUNDAMENTAL);
1732        let expected = 4.0 / 3.0;
1733        assert!(
1734            (c2 - expected).abs() < 1e-10,
1735            "Casimir of fundamental should be 4/3, got {}",
1736            c2
1737        );
1738    }
1739
1740    #[test]
1741    fn test_su3_casimir_antifundamental() {
1742        // (0,1) antifundamental: c₂ = 4/3
1743        use crate::Casimir;
1744        use crate::Su3Irrep;
1745
1746        let c2 = Su3Algebra::quadratic_casimir_eigenvalue(&Su3Irrep::ANTIFUNDAMENTAL);
1747        let expected = 4.0 / 3.0;
1748        assert!(
1749            (c2 - expected).abs() < 1e-10,
1750            "Casimir of antifundamental should be 4/3, got {}",
1751            c2
1752        );
1753    }
1754
1755    #[test]
1756    fn test_su3_casimir_adjoint() {
1757        // (1,1) adjoint (gluon): c₂ = 3
1758        use crate::Casimir;
1759        use crate::Su3Irrep;
1760
1761        let c2 = Su3Algebra::quadratic_casimir_eigenvalue(&Su3Irrep::ADJOINT);
1762        assert_eq!(c2, 3.0, "Casimir of adjoint representation should be 3");
1763    }
1764
1765    #[test]
1766    fn test_su3_casimir_symmetric() {
1767        // (2,0) symmetric (diquark): c₂ = 10/3
1768        use crate::Casimir;
1769        use crate::Su3Irrep;
1770
1771        let c2 = Su3Algebra::quadratic_casimir_eigenvalue(&Su3Irrep::SYMMETRIC);
1772        let expected = 10.0 / 3.0;
1773        assert!(
1774            (c2 - expected).abs() < 1e-10,
1775            "Casimir of symmetric should be 10/3, got {}",
1776            c2
1777        );
1778    }
1779
1780    #[test]
1781    fn test_su3_casimir_formula() {
1782        // Test the formula c₂(p,q) = (1/3)(p² + q² + pq + 3p + 3q)
1783        use crate::Casimir;
1784        use crate::Su3Irrep;
1785
1786        for p in 0..5 {
1787            for q in 0..5 {
1788                let irrep = Su3Irrep::new(p, q);
1789                let c2 = Su3Algebra::quadratic_casimir_eigenvalue(&irrep);
1790
1791                let pf = p as f64;
1792                let qf = q as f64;
1793                let expected = (pf * pf + qf * qf + pf * qf + 3.0 * pf + 3.0 * qf) / 3.0;
1794
1795                assert!(
1796                    (c2 - expected).abs() < 1e-10,
1797                    "Casimir for ({},{}) should be {}, got {}",
1798                    p,
1799                    q,
1800                    expected,
1801                    c2
1802                );
1803            }
1804        }
1805    }
1806
1807    #[test]
1808    fn test_su3_rank() {
1809        // SU(3) has rank 2
1810        use crate::Casimir;
1811
1812        assert_eq!(Su3Algebra::rank(), 2, "SU(3) should have rank 2");
1813        assert_eq!(
1814            Su3Algebra::num_casimirs(),
1815            2,
1816            "SU(3) should have 2 Casimir operators"
1817        );
1818    }
1819
1820    // ==========================================================================
1821    // Cubic Casimir Tests
1822    // ==========================================================================
1823
1824    #[test]
1825    fn test_su3_cubic_casimir_trivial() {
1826        // (0,0) trivial: c₃ = 0
1827        use crate::Casimir;
1828        use crate::Su3Irrep;
1829
1830        let c3_vec = Su3Algebra::higher_casimir_eigenvalues(&Su3Irrep::TRIVIAL);
1831        assert_eq!(c3_vec.len(), 1, "Should return exactly one higher Casimir");
1832        assert_eq!(c3_vec[0], 0.0, "Cubic Casimir of trivial should be 0");
1833    }
1834
1835    #[test]
1836    fn test_su3_cubic_casimir_fundamental() {
1837        // (1,0) fundamental: c₃ = (1/18)(1-0)(2+0+3)(1+0+3) = (1/18)(1)(5)(4) = 20/18 = 10/9
1838        use crate::Casimir;
1839        use crate::Su3Irrep;
1840
1841        let c3_vec = Su3Algebra::higher_casimir_eigenvalues(&Su3Irrep::FUNDAMENTAL);
1842        let c3 = c3_vec[0];
1843        let expected = 10.0 / 9.0;
1844        assert!(
1845            (c3 - expected).abs() < 1e-10,
1846            "Cubic Casimir of fundamental should be 10/9, got {}",
1847            c3
1848        );
1849    }
1850
1851    #[test]
1852    fn test_su3_cubic_casimir_antifundamental() {
1853        // (0,1) antifundamental: c₃ = (1/18)(0-1)(0+1+3)(0+2+3) = (1/18)(-1)(4)(5) = -20/18 = -10/9
1854        use crate::Casimir;
1855        use crate::Su3Irrep;
1856
1857        let c3_vec = Su3Algebra::higher_casimir_eigenvalues(&Su3Irrep::ANTIFUNDAMENTAL);
1858        let c3 = c3_vec[0];
1859        let expected = -10.0 / 9.0;
1860        assert!(
1861            (c3 - expected).abs() < 1e-10,
1862            "Cubic Casimir of antifundamental should be -10/9, got {}",
1863            c3
1864        );
1865    }
1866
1867    #[test]
1868    fn test_su3_cubic_casimir_adjoint() {
1869        // (1,1) adjoint: c₃ = (1/18)(1-1)(2+1+3)(1+2+3) = 0 (self-conjugate)
1870        use crate::Casimir;
1871        use crate::Su3Irrep;
1872
1873        let c3_vec = Su3Algebra::higher_casimir_eigenvalues(&Su3Irrep::ADJOINT);
1874        let c3 = c3_vec[0];
1875        assert!(
1876            c3.abs() < 1e-10,
1877            "Cubic Casimir of adjoint (self-conjugate) should be 0, got {}",
1878            c3
1879        );
1880    }
1881
1882    #[test]
1883    fn test_su3_cubic_casimir_symmetric() {
1884        // (2,0) symmetric: c₃ = (1/18)(2-0)(4+0+3)(2+0+3) = (1/18)(2)(7)(5) = 70/18 = 35/9
1885        use crate::Casimir;
1886        use crate::Su3Irrep;
1887
1888        let c3_vec = Su3Algebra::higher_casimir_eigenvalues(&Su3Irrep::SYMMETRIC);
1889        let c3 = c3_vec[0];
1890        let expected = 70.0 / 18.0;
1891        assert!(
1892            (c3 - expected).abs() < 1e-10,
1893            "Cubic Casimir of symmetric should be 70/18, got {}",
1894            c3
1895        );
1896    }
1897
1898    #[test]
1899    fn test_su3_cubic_casimir_conjugation_symmetry() {
1900        // c₃(p,q) = -c₃(q,p) for conjugate representations
1901        use crate::Casimir;
1902        use crate::Su3Irrep;
1903
1904        for p in 0..5 {
1905            for q in 0..5 {
1906                let irrep_pq = Su3Irrep::new(p, q);
1907                let irrep_qp = Su3Irrep::new(q, p);
1908
1909                let c3_pq = Su3Algebra::higher_casimir_eigenvalues(&irrep_pq)[0];
1910                let c3_qp = Su3Algebra::higher_casimir_eigenvalues(&irrep_qp)[0];
1911
1912                assert!(
1913                    (c3_pq + c3_qp).abs() < 1e-10,
1914                    "c₃({},{}) = {} should equal -c₃({},{}) = {}",
1915                    p,
1916                    q,
1917                    c3_pq,
1918                    q,
1919                    p,
1920                    -c3_qp
1921                );
1922            }
1923        }
1924    }
1925
1926    #[test]
1927    fn test_su3_cubic_casimir_formula() {
1928        // Test the formula c₃(p,q) = (1/18)(p-q)(2p+q+3)(p+2q+3)
1929        use crate::Casimir;
1930        use crate::Su3Irrep;
1931
1932        for p in 0..5 {
1933            for q in 0..5 {
1934                let irrep = Su3Irrep::new(p, q);
1935                let c3 = Su3Algebra::higher_casimir_eigenvalues(&irrep)[0];
1936
1937                let pf = p as f64;
1938                let qf = q as f64;
1939                let expected = (pf - qf) * (2.0 * pf + qf + 3.0) * (pf + 2.0 * qf + 3.0) / 18.0;
1940
1941                assert!(
1942                    (c3 - expected).abs() < 1e-10,
1943                    "Cubic Casimir for ({},{}) should be {}, got {}",
1944                    p,
1945                    q,
1946                    expected,
1947                    c3
1948                );
1949            }
1950        }
1951    }
1952
1953    // ==========================================================================
1954    // Edge Case Tests: Near-Singular Scenarios
1955    // ==========================================================================
1956    //
1957    // These tests verify numerical stability at boundary conditions, particularly
1958    // for the scale-relative threshold improvements.
1959
1960    /// Test Gram-Schmidt with very small matrices
1961    ///
1962    /// Ensures the scale-relative threshold handles matrices with small norms
1963    /// without false positives for numerical breakdown.
1964    #[test]
1965    fn test_gram_schmidt_small_matrix() {
1966        // Create a matrix with small norm (scale ~ 1e-6)
1967        let small_algebra = Su3Algebra([1e-6, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]);
1968        let small_element = SU3::exp(&small_algebra);
1969
1970        // Should still be unitary
1971        assert!(
1972            small_element.verify_unitarity(1e-10),
1973            "Small element should be unitary"
1974        );
1975
1976        // exp followed by log should round-trip
1977        let recovered = small_element
1978            .log()
1979            .expect("log should succeed for near-identity");
1980        let diff = (recovered.0[0] - small_algebra.0[0]).abs();
1981        assert!(
1982            diff < 1e-10,
1983            "exp/log round-trip failed for small element: diff = {}",
1984            diff
1985        );
1986    }
1987
1988    /// Test Gram-Schmidt with very large matrices
1989    ///
1990    /// Ensures the scale-relative threshold handles matrices with large norms
1991    /// correctly, where absolute thresholds would be too tight.
1992    #[test]
1993    fn test_gram_schmidt_large_matrix() {
1994        use crate::traits::LieGroup;
1995
1996        // Create a matrix with large rotation angle (near π)
1997        let large_algebra = Su3Algebra([2.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]);
1998        let large_element = SU3::exp(&large_algebra);
1999
2000        // Should still be unitary after exp (reorthogonalization should help)
2001        assert!(
2002            large_element.verify_unitarity(1e-8),
2003            "Large rotation element should be unitary, distance = {}",
2004            large_element.distance_to_identity()
2005        );
2006    }
2007
2008    /// Test exp with repeated squaring stability
2009    ///
2010    /// Verifies that reorthogonalization every squaring prevents drift.
2011    #[test]
2012    fn test_exp_repeated_squaring_stability() {
2013        use crate::traits::LieGroup;
2014
2015        // Angle that requires multiple squaring steps (2^k scaling)
2016        // θ = 3.0 requires about 3-4 squaring steps
2017        let algebra = Su3Algebra([1.0, 1.0, 0.5, 0.3, 0.2, 0.1, 0.0, 0.0]);
2018        let element = SU3::exp(&algebra);
2019
2020        // Check unitarity is preserved
2021        assert!(
2022            element.verify_unitarity(1e-10),
2023            "Repeated squaring should preserve unitarity"
2024        );
2025
2026        // Check special unitarity: U†U = I implies det(U) has magnitude 1
2027        // For SU(3), we also require det = 1 (not just |det| = 1)
2028        // We verify this indirectly through verify_special_unitarity if available,
2029        // or by checking trace(U†U) = 3 (trace of identity)
2030        let u_dag_u = element.adjoint().compose(&element);
2031        let trace = u_dag_u.trace();
2032        assert!(
2033            (trace.re - 3.0).abs() < 1e-10 && trace.im.abs() < 1e-10,
2034            "U†U should have trace 3, got {:.6}+{:.6}i",
2035            trace.re,
2036            trace.im
2037        );
2038    }
2039
2040    /// Test exp/log round-trip near group boundary
2041    ///
2042    /// For elements near θ = π, log can have numerical difficulties.
2043    /// This test verifies graceful handling.
2044    #[test]
2045    fn test_exp_log_near_boundary() {
2046        use crate::traits::LieGroup;
2047
2048        // Angle close to but not exactly at π
2049        let theta = std::f64::consts::PI - 0.1;
2050        let algebra = Su3Algebra([theta, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]);
2051        let element = SU3::exp(&algebra);
2052
2053        // Element should be valid
2054        assert!(
2055            element.verify_unitarity(1e-10),
2056            "Near-boundary element should be unitary"
2057        );
2058
2059        // Log might fail or give approximate result - just verify it doesn't panic
2060        if let Ok(recovered) = element.log() {
2061            // If log succeeds, verify exp(log(g)) = g
2062            let round_trip = SU3::exp(&recovered);
2063            let distance = element.distance(&round_trip);
2064            assert!(
2065                distance < 1e-6,
2066                "Round trip should be close, distance = {}",
2067                distance
2068            );
2069        }
2070        // Log failure near θ=π is acceptable
2071        // The element is still valid, just log has a singularity
2072    }
2073
2074    /// Test composition preserves unitarity under repeated operations
2075    ///
2076    /// Accumulated errors from many compositions could drift from SU(3).
2077    #[test]
2078    fn test_composition_stability() {
2079        use crate::traits::LieGroup;
2080
2081        // Create several distinct SU(3) elements
2082        let u1 = SU3::exp(&Su3Algebra([0.5, 0.3, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]));
2083        let u2 = SU3::exp(&Su3Algebra([0.0, 0.4, 0.2, 0.0, 0.0, 0.0, 0.0, 0.0]));
2084        let u3 = SU3::exp(&Su3Algebra([0.0, 0.0, 0.3, 0.5, 0.0, 0.0, 0.0, 0.0]));
2085
2086        // Compose many times
2087        let mut result = SU3::identity();
2088        for _ in 0..100 {
2089            result = result.compose(&u1);
2090            result = result.compose(&u2);
2091            result = result.compose(&u3);
2092        }
2093
2094        // Should still be unitary (SU(3) is closed under composition)
2095        assert!(
2096            result.verify_unitarity(1e-8),
2097            "100 compositions should still be unitary"
2098        );
2099    }
2100}