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