Skip to main content

lie_groups/
traits.rs

1//! Traits for Lie groups and Lie algebras.
2//!
3//! This module defines the core abstractions for Lie theory. These traits
4//! capture the mathematical structure
5//! of Lie groups as both algebraic objects (groups) and geometric objects
6//! (differentiable manifolds).
7//!
8//! # Mathematical Background
9//!
10//! A **Lie group** G is a smooth manifold that is also a group, where the
11//! group operations are smooth maps:
12//!
13//! - Multiplication: G × G → G is smooth
14//! - Inversion: G → G is smooth
15//! - Identity: e ∈ G
16//!
17//! A **Lie algebra** 𝔤 is the tangent space at the identity, with a Lie bracket
18//! operation [·,·] : 𝔤 × 𝔤 → 𝔤 satisfying bilinearity, antisymmetry, and the Jacobi identity.
19//!
20//! The **exponential map** exp: 𝔤 → G connects the Lie algebra to the Lie group,
21//! mapping tangent vectors to group elements via one-parameter subgroups.
22//!
23//! # Design Philosophy
24//!
25//! These traits enable **generic programming over Lie groups**: algorithms written
26//! against the `LieGroup` trait work for SU(2), U(1), SU(3), etc. without modification.
27//! This mirrors how mathematicians work—proving theorems for "a Lie group G"
28//! rather than for each specific group.
29//!
30//! # Examples
31//!
32//! ```
33//! use lie_groups::{LieGroup, SU2};
34//!
35//! // Generic parallel transport along a path
36//! fn parallel_transport<G: LieGroup>(path: &[G]) -> G {
37//!     path.iter().fold(G::identity(), |acc, g| acc.compose(g))
38//! }
39//!
40//! // Works for any Lie group!
41//! let su2_path = vec![SU2::rotation_x(0.1), SU2::rotation_y(0.2)];
42//! let holonomy = parallel_transport(&su2_path);
43//! ```
44//!
45//! # Constrained Representations: Compile-Time Guarantees
46//!
47//! This module uses **constrained representations** to encode mathematical
48//! properties at the type level. The key insight: instead of checking if a
49//! value satisfies a property at runtime, use a representation that makes
50//! invalid values **unrepresentable**.
51//!
52//! ## The Pattern
53//!
54//! | Property | Representation | Why It Works |
55//! |----------|----------------|--------------|
56//! | Traceless | Basis coefficients | Pauli/Gell-Mann matrices are traceless |
57//! | Anti-Hermitian | `i·(Hermitian basis)` | `(iH)† = -iH† = -iH` |
58//! | Unit norm | Private fields + normalizing constructor | Can only create unit values |
59//! | Fixed size | `[T; N]` instead of `Vec<T>` | Compile-time size guarantee |
60//!
61//! ## Example: Traceless by Construction
62//!
63//! ```text
64//! // BAD: General matrix, might not be traceless
65//! struct BadAlgebra { matrix: [[f64; 2]; 2] }  // tr(M) could be anything
66//!
67//! // GOOD: Pauli basis coefficients, always traceless
68//! struct Su2Algebra([f64; 3]);  // X = a·iσ₁ + b·iσ₂ + c·iσ₃, tr(X) = 0 always
69//! ```
70//!
71//! ## Connection to Lean Proofs
72//!
73//! The marker traits in this module connect to formal proofs:
74//!
75//! - [`TracelessByConstruction`]: Combined with the Lean theorem
76//!   `det_exp_eq_exp_trace`, gives `det(exp(X)) = exp(0) = 1` for all X.
77//!
78//! - [`AntiHermitianByConstruction`]: Combined with `exp_antiHermitian_unitary`,
79//!   gives `exp(X)† · exp(X) = I` for all X.
80//!
81//! These are **compile-time guarantees**: if your code compiles with the trait
82//! bound `A: TracelessByConstruction`, you know all algebra elements are
83//! traceless without any runtime checks.
84//!
85//! ## Zero Runtime Overhead
86//!
87//! Marker traits are empty—no methods, no data:
88//!
89//! ```ignore
90//! pub trait TracelessByConstruction: LieAlgebra {}
91//! impl TracelessByConstruction for Su2Algebra {}  // Empty impl
92//! ```
93//!
94//! The compiler uses them for type checking, then erases them completely.
95//! Same pattern as `Send`, `Sync`, `Copy` in std.
96//!
97//! ## When to Use Each Pattern
98//!
99//! | Pattern | Use When | Example |
100//! |---------|----------|---------|
101//! | Sealed marker trait | Property depends on type, not value | `TracelessByConstruction` |
102//! | Private fields | Property depends on value (invariant) | `UnitQuaternion` unit norm |
103//! | Fixed-size array | Size is known at compile time | `Plaquette([NodeIndex; 4])` |
104//! | Builder pattern | Invariant built incrementally | Path with consecutive edges |
105
106// ============================================================================
107// Sealed Trait Pattern for Marker Traits
108// ============================================================================
109//
110// Marker traits (Compact, Abelian, Simple, SemiSimple) are sealed to prevent
111// incorrect claims about mathematical properties. These encode theorems:
112// - "U(1) is abelian" is a mathematical fact that shouldn't be claimable by arbitrary types
113// - "SU(2) is compact" likewise
114//
115// The core `LieGroup` and `LieAlgebra` traits are OPEN - the associated type
116// constraints provide compile-time safety without sealing.
117
118mod sealed {
119    // Import types for cleaner impl blocks
120    use super::super::rplus::RPlus;
121    use super::super::so3::{So3Algebra, SO3};
122    use super::super::su2::{Su2Algebra, SU2};
123    use super::super::su3::{Su3Algebra, SU3};
124    use super::super::sun::{SunAlgebra, SUN};
125    use super::super::u1::{U1Algebra, U1};
126
127    /// Sealed trait for compact groups.
128    pub trait SealedCompact {}
129
130    /// Sealed trait for abelian groups.
131    pub trait SealedAbelian {}
132
133    /// Sealed trait for simple groups.
134    pub trait SealedSimple {}
135
136    /// Sealed trait for semi-simple groups.
137    pub trait SealedSemiSimple {}
138
139    /// Sealed trait for algebras that are traceless by construction.
140    ///
141    /// An algebra is "traceless by construction" if its representation
142    /// structurally guarantees tr(X) = 0 for all elements X.
143    ///
144    /// # Examples
145    /// - `Su2Algebra::new([f64; 3])`: Pauli basis coefficients → always traceless
146    /// - `Su3Algebra::new([f64; 8])`: Gell-Mann basis coefficients → always traceless
147    ///
148    /// # Non-examples
149    /// - `U1Algebra::new(f64)`: Represents iℝ, trace = i·a ≠ 0
150    /// - General matrix algebras: Can have arbitrary trace
151    pub trait SealedTraceless {}
152
153    /// Sealed trait for algebras that are anti-Hermitian by construction.
154    ///
155    /// An algebra is "anti-Hermitian by construction" if its representation
156    /// structurally guarantees X† = -X for all elements X.
157    ///
158    /// # Examples
159    /// - `Su2Algebra`: Uses i·σ basis (anti-Hermitian)
160    /// - `U1Algebra`: Represents iℝ ⊂ ℂ (purely imaginary = anti-Hermitian)
161    /// - `So3Algebra`: Real antisymmetric matrices (= anti-Hermitian)
162    pub trait SealedAntiHermitian {}
163
164    // =========================================================================
165    // Sealed trait implementations
166    // =========================================================================
167
168    // Compact groups (R⁺ is non-compact, hence excluded)
169    impl SealedCompact for SU2 {}
170    impl SealedCompact for SU3 {}
171    impl SealedCompact for U1 {}
172    impl SealedCompact for SO3 {}
173    impl<const N: usize> SealedCompact for SUN<N> {}
174
175    // Abelian groups
176    impl SealedAbelian for U1 {}
177    impl SealedAbelian for RPlus {}
178
179    // Simple groups (SU(2), SU(3), SO(3), SU(N) for N >= 2)
180    impl SealedSimple for SU2 {}
181    impl SealedSimple for SU3 {}
182    impl SealedSimple for SO3 {}
183    // SUN<N> is simple for N >= 2. SU(1) is trivial (not simple), but Rust
184    // stable lacks const generic bounds (`where N >= 2`). The runtime assert
185    // in SUN::identity/exp/from_matrix prevents SUN<1> construction, so the
186    // marker is unreachable in practice. Revisit when `generic_const_exprs`
187    // stabilizes.
188    impl<const N: usize> SealedSimple for SUN<N> {}
189
190    // Semi-simple groups (all simple groups are semi-simple)
191    impl SealedSemiSimple for SU2 {}
192    impl SealedSemiSimple for SU3 {}
193    impl SealedSemiSimple for SO3 {}
194    // Same N >= 2 caveat as SealedSimple above.
195    impl<const N: usize> SealedSemiSimple for SUN<N> {}
196
197    // Traceless algebras (representation guarantees tr(X) = 0)
198    // NOT U1Algebra: u(1) = iℝ has tr(ia) = ia ≠ 0
199    impl SealedTraceless for Su2Algebra {}
200    impl SealedTraceless for Su3Algebra {}
201    impl SealedTraceless for So3Algebra {}
202    impl<const N: usize> SealedTraceless for SunAlgebra<N> {}
203
204    // Anti-Hermitian algebras (representation guarantees X† = -X)
205    impl SealedAntiHermitian for Su2Algebra {}
206    impl SealedAntiHermitian for Su3Algebra {}
207    impl SealedAntiHermitian for So3Algebra {}
208    impl SealedAntiHermitian for U1Algebra {}
209    impl<const N: usize> SealedAntiHermitian for SunAlgebra<N> {}
210}
211
212// ============================================================================
213// Marker Traits for Mathematical Properties
214// ============================================================================
215
216/// Marker trait for compact Lie groups.
217///
218/// A Lie group is **compact** if it is compact as a topological space.
219///
220/// # Examples
221/// - `U(1)`, `SU(2)`, `SU(3)`, `SO(3)` - all compact
222///
223/// # Significance
224/// - Representations are completely reducible
225/// - Yang-Mills theory well-defined
226/// - Haar measure exists
227///
228/// # Sealed Trait
229///
230/// This trait is sealed - only verified compact groups can implement it.
231pub trait Compact: LieGroup + sealed::SealedCompact {}
232
233/// Marker trait for abelian (commutative) Lie groups.
234///
235/// Satisfies: `g.compose(&h) == h.compose(&g)` for all g, h
236///
237/// # Examples
238/// - `U(1)` - abelian (electromagnetism)
239///
240/// # Type Safety
241/// ```ignore
242/// fn maxwell_theory<G: Abelian>(conn: &NetworkConnection<G>) {
243///     // Compiles only for abelian groups!
244/// }
245/// ```
246///
247/// # Sealed Trait
248///
249/// This trait is sealed - only verified abelian groups can implement it.
250pub trait Abelian: LieGroup + sealed::SealedAbelian {}
251
252/// Marker trait for simple Lie groups (no non-trivial normal subgroups).
253///
254/// # Examples
255/// - `SU(2)`, `SU(3)`, `SO(3)` - simple
256///
257/// # Sealed Trait
258///
259/// This trait is sealed - only verified simple groups can implement it.
260pub trait Simple: SemiSimple + sealed::SealedSimple {}
261
262/// Marker trait for semi-simple Lie groups (products of simple groups).
263///
264/// # Examples
265/// - `SU(2)`, `SU(3)` - semi-simple
266///
267/// # Sealed Trait
268///
269/// This trait is sealed - only verified semi-simple groups can implement it.
270pub trait SemiSimple: LieGroup + sealed::SealedSemiSimple {}
271
272/// Marker trait for Lie algebras whose representation is structurally traceless.
273///
274/// A Lie algebra is "traceless by construction" when its representation
275/// **guarantees** tr(X) = 0 for all elements X, without runtime verification.
276///
277/// # Mathematical Significance
278///
279/// Combined with the theorem `det(exp(X)) = exp(tr(X))` (proven in Lean as
280/// `Matrix.det_exp_eq_exp_trace`), this gives a **compile-time guarantee**:
281///
282/// ```text
283/// A: TracelessByConstruction  ⟹  tr(X) = 0 for all X ∈ A
284///                             ⟹  det(exp(X)) = exp(0) = 1
285/// ```
286///
287/// Therefore, `exp` maps traceless algebras to special (det = 1) groups.
288///
289/// # Examples
290///
291/// - `Su2Algebra::new([f64; 3])`: Pauli basis coefficients → always traceless
292/// - `Su3Algebra::new([f64; 8])`: Gell-Mann basis coefficients → always traceless
293/// - `So3Algebra::new([f64; 3])`: Antisymmetric matrices → always traceless
294///
295/// # Non-examples
296///
297/// - `U1Algebra::new(f64)`: Represents iℝ ⊂ ℂ, where tr(ia) = ia ≠ 0
298/// - General `Matrix<N,N>`: Can have arbitrary trace
299///
300/// # Type-Level Proof Pattern
301///
302/// ```ignore
303/// fn exp_to_special_unitary<A>(x: &A) -> SpecialUnitary
304/// where
305///     A: TracelessByConstruction,
306/// {
307///     // Compiler knows: tr(x) = 0 by construction
308///     // Lean theorem: det(exp(x)) = exp(tr(x)) = exp(0) = 1
309///     // Therefore: exp(x) ∈ SU(n), not just U(n)
310/// }
311/// ```
312///
313/// # Sealed Trait
314///
315/// This trait is sealed - only algebras with verified traceless representations
316/// can implement it.
317pub trait TracelessByConstruction: LieAlgebra + sealed::SealedTraceless {}
318
319/// Marker trait for Lie algebras whose representation is structurally anti-Hermitian.
320///
321/// A Lie algebra is "anti-Hermitian by construction" when its representation
322/// **guarantees** X† = -X for all elements X, without runtime verification.
323///
324/// # Mathematical Significance
325///
326/// Combined with the theorem `exp(X)† · exp(X) = I` for anti-Hermitian X
327/// (proven in Lean as `Matrix.exp_antiHermitian_unitary`), this gives a
328/// **compile-time guarantee**:
329///
330/// ```text
331/// A: AntiHermitianByConstruction  ⟹  X† = -X for all X ∈ A
332///                                 ⟹  exp(X) is unitary
333/// ```
334///
335/// # Examples
336///
337/// - `Su2Algebra`: Uses i·σ basis (anti-Hermitian Pauli matrices)
338/// - `Su3Algebra`: Uses i·λ basis (anti-Hermitian Gell-Mann matrices)
339/// - `So3Algebra`: Real antisymmetric matrices (= anti-Hermitian over ℝ)
340/// - `U1Algebra`: Represents iℝ ⊂ ℂ (purely imaginary = anti-Hermitian)
341///
342/// # Sealed Trait
343///
344/// This trait is sealed - only algebras with verified anti-Hermitian representations
345/// can implement it.
346pub trait AntiHermitianByConstruction: LieAlgebra + sealed::SealedAntiHermitian {}
347
348/// Trait for Lie algebra elements.
349///
350/// A Lie algebra 𝔤 is a vector space equipped with a bilinear, antisymmetric
351/// operation called the Lie bracket [·,·] : 𝔤 × 𝔤 → 𝔤 that satisfies the Jacobi identity.
352///
353/// In this library, we represent Lie algebras as the tangent space at the identity
354/// of a Lie group. For matrix groups, algebra elements are traceless anti-Hermitian matrices.
355///
356/// # Mathematical Structure
357///
358/// For classical Lie groups:
359/// - **u(1)**: ℝ (1-dimensional)
360/// - **su(2)**: ℝ³ (3-dimensional, isomorphic to pure quaternions)
361/// - **su(3)**: ℝ⁸ (8-dimensional, Gell-Mann basis)
362///
363/// # Design
364///
365/// This trait enables generic gradient descent on Lie groups via the exponential map.
366/// The gradient lives in the Lie algebra, and we use `exp: 𝔤 → G` to move on the manifold.
367///
368/// # Examples
369///
370/// ```ignore
371/// use lie_groups::{LieAlgebra, Su2Algebra};
372///
373/// let v = Su2Algebra::basis_element(0);  // Pauli X direction
374/// let w = Su2Algebra::basis_element(1);  // Pauli Y direction
375/// let sum = v.add(&w);                   // Linear combination
376/// let scaled = v.scale(2.0);             // Scalar multiplication
377/// ```
378///
379/// # Open Trait
380///
381/// This trait is open for implementation. The associated type constraints
382/// on `LieGroup` (requiring `type Algebra: LieAlgebra`) provide compile-time
383/// safety by enforcing group-algebra correspondence.
384pub trait LieAlgebra: Clone + Sized + std::fmt::Debug + PartialEq {
385    /// Dimension of the Lie algebra as a compile-time constant.
386    ///
387    /// This enables compile-time dimension checking where possible.
388    /// For example, `from_components` can verify array length at compile time
389    /// when called with a fixed-size array.
390    ///
391    /// # Values
392    ///
393    /// - `U1Algebra`: 1
394    /// - `Su2Algebra`: 3
395    /// - `So3Algebra`: 3
396    /// - `Su3Algebra`: 8
397    /// - `SunAlgebra<N>`: N² - 1
398    ///
399    /// # Note
400    ///
401    /// Rust doesn't yet support `where A: LieAlgebra<DIM = 3>` bounds
402    /// (requires `generic_const_exprs`). Use `const_assert!` or runtime
403    /// checks with `Self::DIM` for now.
404    const DIM: usize;
405
406    /// Zero element (additive identity) 0 ∈ 𝔤.
407    ///
408    /// Satisfies: `v.add(&Self::zero()) == v` for all v ∈ 𝔤.
409    #[must_use]
410    fn zero() -> Self;
411
412    /// Add two algebra elements: v + w
413    ///
414    /// The Lie algebra is a vector space, so addition is defined.
415    ///
416    /// # Properties
417    ///
418    /// - Commutative: `v.add(&w) == w.add(&v)`
419    /// - Associative: `u.add(&v.add(&w)) == u.add(&v).add(&w)`
420    /// - Identity: `v.add(&Self::zero()) == v`
421    #[must_use]
422    fn add(&self, other: &Self) -> Self;
423
424    /// Scalar multiplication: α · v
425    ///
426    /// Scale the algebra element by a real number.
427    ///
428    /// # Properties
429    ///
430    /// - Distributive: `v.scale(a + b) ≈ v.scale(a).add(&v.scale(b))`
431    /// - Associative: `v.scale(a * b) ≈ v.scale(a).scale(b)`
432    /// - Identity: `v.scale(1.0) == v`
433    #[must_use]
434    fn scale(&self, scalar: f64) -> Self;
435
436    /// Euclidean norm of the coefficient vector: ||v|| = √(Σᵢ vᵢ²)
437    ///
438    /// This is the L² norm in the **coefficient space**, not the matrix operator
439    /// norm. For the standard Gell-Mann normalization `Tr(λᵢ λⱼ) = 2δᵢⱼ`,
440    /// the coefficient norm equals `||X||_F / √2` where `||·||_F` is the
441    /// Frobenius norm of the matrix representation.
442    ///
443    /// # Convergence Bounds
444    ///
445    /// The BCH convergence radius `||X||_op + ||Y||_op < log(2)` uses the
446    /// operator norm. Since `||X||_op ≤ ||X||_F` and `||X||_F` is proportional
447    /// to the coefficient norm (with a constant depending on basis normalization),
448    /// `x.norm() + y.norm() < log(2)` provides a sufficient (conservative)
449    /// convergence check for standard normalizations (tr(TₐTᵦ) = ½δₐᵦ).
450    ///
451    /// # Returns
452    ///
453    /// Non-negative real number. Zero iff `self` is the zero element.
454    #[must_use]
455    fn norm(&self) -> f64;
456
457    /// Get the i-th basis element of the Lie algebra.
458    ///
459    /// For su(2), the basis is {σₓ, σᵧ, σᵤ}/2i (Pauli matrices).
460    /// For u(1), the single basis element is i (imaginary unit).
461    ///
462    /// # Parameters
463    ///
464    /// - `i`: Index in [0, `dim()`)
465    ///
466    /// # Panics
467    ///
468    /// If `i >= Self::DIM`.
469    ///
470    /// # Example
471    ///
472    /// ```ignore
473    /// use lie_groups::{LieAlgebra, Su2Algebra};
474    ///
475    /// let sigma_x = Su2Algebra::basis_element(0);
476    /// let sigma_y = Su2Algebra::basis_element(1);
477    /// let sigma_z = Su2Algebra::basis_element(2);
478    /// ```
479    #[must_use]
480    fn basis_element(i: usize) -> Self;
481
482    /// Construct algebra element from basis coordinates.
483    ///
484    /// Given coefficients [c₀, c₁, ..., c_{n-1}], constructs:
485    /// ```text
486    /// v = Σᵢ cᵢ eᵢ
487    /// ```
488    /// where eᵢ are the basis elements.
489    ///
490    /// # Parameters
491    ///
492    /// - `components`: Slice of length `Self::DIM`
493    ///
494    /// # Panics
495    ///
496    /// If `components.len() != Self::DIM`.
497    ///
498    /// # Example
499    ///
500    /// ```ignore
501    /// use lie_groups::{LieAlgebra, Su2Algebra};
502    ///
503    /// let v = Su2Algebra::from_components(&[1.0, 0.0, 0.0]);  // Pure X
504    /// let w = Su2Algebra::from_components(&[0.5, 0.5, 0.0]);  // X + Y mix
505    /// ```
506    #[must_use]
507    fn from_components(components: &[f64]) -> Self;
508
509    /// Extract basis coordinates from algebra element.
510    ///
511    /// Returns the coefficients [c₀, c₁, ..., c_{n-1}] such that:
512    /// ```text
513    /// self = Σᵢ cᵢ eᵢ
514    /// ```
515    /// where eᵢ are the basis elements.
516    ///
517    /// # Returns
518    ///
519    /// Vector of length `Self::DIM` containing the basis coordinates.
520    ///
521    /// # Example
522    ///
523    /// ```ignore
524    /// use lie_groups::{LieAlgebra, Su2Algebra};
525    ///
526    /// let v = Su2Algebra::from_components(&[1.0, 2.0, 3.0]);
527    /// let coords = v.to_components();
528    /// assert_eq!(coords, vec![1.0, 2.0, 3.0]);
529    /// ```
530    #[must_use]
531    fn to_components(&self) -> Vec<f64>;
532
533    /// Lie bracket operation: [X, Y] ∈ 𝔤
534    ///
535    /// The Lie bracket is a binary operation on the Lie algebra that captures
536    /// the non-commutativity of the group. For matrix Lie algebras:
537    /// ```text
538    /// [X, Y] = XY - YX  (matrix commutator)
539    /// ```
540    ///
541    /// # Mathematical Properties
542    ///
543    /// The Lie bracket must satisfy:
544    ///
545    /// 1. **Bilinearity**:
546    ///    - `[aX + bY, Z] = a[X,Z] + b[Y,Z]`
547    ///    - `[X, aY + bZ] = a[X,Y] + b[X,Z]`
548    ///
549    /// 2. **Antisymmetry**: `[X, Y] = -[Y, X]`
550    ///
551    /// 3. **Jacobi Identity**:
552    ///    ```text
553    ///    [X, [Y, Z]] + [Y, [Z, X]] + [Z, [X, Y]] = 0
554    ///    ```
555    ///
556    /// # Geometric Interpretation
557    ///
558    /// The bracket measures the "infinitesimal non-commutativity" of the group.
559    /// For SU(2) with structure constants `[eᵢ, eⱼ] = −εᵢⱼₖeₖ`:
560    /// ```text
561    /// [e₁, e₂] = −e₃,  [e₂, e₃] = −e₁,  [e₃, e₁] = −e₂
562    /// ```
563    ///
564    /// # Applications
565    ///
566    /// - **Structure constants**: Define the algebra's multiplication table
567    /// - **Curvature formulas**: F = dA + [A, A]
568    /// - **Representation theory**: Adjoint representation
569    ///
570    /// # Examples
571    ///
572    /// ```ignore
573    /// use lie_groups::{LieAlgebra, Su2Algebra};
574    ///
575    /// let e1 = Su2Algebra::basis_element(0);
576    /// let e2 = Su2Algebra::basis_element(1);
577    /// let bracket = e1.bracket(&e2);
578    ///
579    /// // [e₁, e₂] = −e₃, so bracket has norm 1
580    /// assert!((bracket.norm() - 1.0).abs() < 1e-10);
581    /// ```
582    ///
583    #[must_use]
584    fn bracket(&self, other: &Self) -> Self;
585
586    /// Inner product on coefficient space: ⟨v, w⟩ = Σᵢ vᵢ wᵢ
587    ///
588    /// This is the Euclidean dot product on basis coefficients, **not** the
589    /// Killing form B(X,Y) = tr(ad X · ad Y). The two are proportional for
590    /// semisimple algebras (e.g., B = −4⟨·,·⟩ for su(2) with the iσ/2 basis),
591    /// but this method always returns the coefficient dot product.
592    ///
593    /// # Default Implementation
594    ///
595    /// Computes `Σᵢ self.to_components()[i] * other.to_components()[i]`.
596    /// Override for algebras with more efficient direct computation.
597    #[must_use]
598    fn inner(&self, other: &Self) -> f64 {
599        let v = self.to_components();
600        let w = other.to_components();
601        v.iter().zip(w.iter()).map(|(vi, wi)| vi * wi).sum()
602    }
603}
604
605/// Core trait for Lie group elements.
606///
607/// A Lie group is a smooth manifold with a compatible group structure.
608/// This trait captures both the algebraic (group operations) and geometric
609/// (distance, smoothness) aspects.
610///
611/// # Mathematical Properties
612///
613/// Implementations must satisfy the group axioms:
614///
615/// 1. **Identity**: `g.compose(&G::identity()) == g`
616/// 2. **Inverse**: `g.compose(&g.inverse()) == G::identity()`
617/// 3. **Associativity**: `g1.compose(&g2.compose(&g3)) == g1.compose(&g2).compose(&g3)`
618/// 4. **Closure**: `compose` always produces a valid group element
619///
620/// Additionally, the manifold structure should satisfy:
621///
622/// 5. **Metric**: `distance_to_identity()` measures geodesic distance on the manifold
623/// 6. **Smoothness**: Small parameter changes produce smooth paths in the group
624///
625/// # Type Parameters vs Associated Types
626///
627/// This trait uses `Self` for the group element type, allowing different
628/// representations of the same abstract group (e.g., matrix vs quaternion
629/// representation of SU(2)).
630///
631/// # Examples
632///
633/// ```
634/// use lie_groups::{LieGroup, SU2};
635///
636/// let g = SU2::rotation_x(1.0);
637/// let h = SU2::rotation_y(0.5);
638///
639/// // Group operations
640/// let product = g.compose(&h);
641/// let inv = g.inverse();
642///
643/// // Check group axioms
644/// assert!(g.compose(&g.inverse()).distance_to_identity() < 1e-10);
645/// ```
646///
647/// # Design Notes
648///
649/// - **`conjugate_transpose` assumes a matrix representation** with an inner
650///   product (unitary/orthogonal groups). For non-unitary groups this method
651///   would need a different interpretation. All groups in this crate are
652///   unitary or orthogonal, so this is not currently a problem.
653///
654/// - **`log` returns `LogResult`** because the logarithm has group-specific
655///   singularities (e.g., -I for SU(2), π-rotation for SO(3)). Some groups
656///   provide richer inherent methods like `log_with_condition()` for
657///   numerical conditioning information.
658///
659/// # Open Trait
660///
661/// This trait is open for implementation. The associated type `Algebra`
662/// must implement `LieAlgebra`, providing compile-time enforcement of
663/// the group-algebra correspondence.
664pub trait LieGroup: Clone + Sized + std::fmt::Debug {
665    /// Matrix dimension in the fundamental representation.
666    ///
667    /// For matrix Lie groups, this is the size of the N×N matrices:
668    /// - U(1): `MATRIX_DIM = 1` (complex numbers as 1×1 matrices)
669    /// - SU(2): `MATRIX_DIM = 2` (2×2 unitary matrices)
670    /// - SU(3): `MATRIX_DIM = 3` (3×3 unitary matrices)
671    /// - SO(3): `MATRIX_DIM = 3` (3×3 orthogonal matrices)
672    ///
673    /// # Distinction from `LieAlgebra::DIM`
674    ///
675    /// `MATRIX_DIM` is the matrix size N, while `LieAlgebra::DIM` is the
676    /// algebra dimension (number of independent generators):
677    /// - SU(2): `MATRIX_DIM = 2`, `Algebra::DIM = 3`
678    /// - SU(3): `MATRIX_DIM = 3`, `Algebra::DIM = 8`
679    const MATRIX_DIM: usize;
680
681    /// Associated Lie algebra type.
682    ///
683    /// The Lie algebra 𝔤 is the tangent space at the identity of the Lie group G.
684    /// This associated type enables generic algorithms that work in the algebra
685    /// (e.g., gradient descent, exponential coordinates).
686    ///
687    /// # Examples
688    ///
689    /// - `SU(2)::Algebra` = `Su2Algebra` (3-dimensional)
690    /// - `U(1)::Algebra` = `U1Algebra` (1-dimensional)
691    /// - `SU(3)::Algebra` = `Su3Algebra` (8-dimensional)
692    type Algebra: LieAlgebra;
693
694    /// The identity element e ∈ G.
695    ///
696    /// Satisfies: `g.compose(&Self::identity()) == g` for all g.
697    ///
698    /// # Mathematical Note
699    ///
700    /// For matrix groups, this is the identity matrix I.
701    /// For additive groups, this would be zero.
702    #[must_use]
703    fn identity() -> Self;
704
705    /// Group composition (multiplication): g₁ · g₂
706    ///
707    /// Computes the product of two group elements. For matrix groups,
708    /// this is matrix multiplication. For SU(2), preserves unitarity.
709    ///
710    /// # Properties
711    ///
712    /// - Associative: `g1.compose(&g2.compose(&g3)) == g1.compose(&g2).compose(&g3)`
713    /// - Identity: `g.compose(&Self::identity()) == g`
714    ///
715    /// # Complexity
716    ///
717    /// | Group | Time | Space | Notes |
718    /// |-------|------|-------|-------|
719    /// | SU(2) via quaternion | O(1) | O(1) | 16 multiplications |
720    /// | SU(N) via matrix | O(N³) | O(N²) | Matrix multiplication |
721    /// | U(1) | O(1) | O(1) | Complex multiplication |
722    ///
723    /// # Examples
724    ///
725    /// ```
726    /// use lie_groups::{LieGroup, SU2};
727    ///
728    /// let g = SU2::rotation_x(0.5);
729    /// let h = SU2::rotation_y(0.3);
730    /// let product = g.compose(&h);
731    /// ```
732    #[must_use]
733    fn compose(&self, other: &Self) -> Self;
734
735    /// Group inverse: g⁻¹
736    ///
737    /// Computes the unique element such that `g.compose(&g.inverse()) == identity`.
738    /// For unitary matrix groups, this is the conjugate transpose (adjoint).
739    ///
740    /// # Properties
741    ///
742    /// - Involutive: `g.inverse().inverse() == g`
743    /// - Inverse property: `g.compose(&g.inverse()) == Self::identity()`
744    ///
745    /// # Complexity
746    ///
747    /// | Group | Time | Space | Notes |
748    /// |-------|------|-------|-------|
749    /// | SU(2) via quaternion | O(1) | O(1) | Quaternion conjugate |
750    /// | SU(N) unitary | O(N²) | O(N²) | Conjugate transpose |
751    /// | U(1) | O(1) | O(1) | Complex conjugate |
752    ///
753    /// # Examples
754    ///
755    /// ```
756    /// use lie_groups::{LieGroup, SU2};
757    ///
758    /// let g = SU2::rotation_z(1.2);
759    /// let g_inv = g.inverse();
760    ///
761    /// // Verify inverse property (up to numerical precision)
762    /// assert!(g.compose(&g_inv).distance_to_identity() < 1e-10);
763    /// ```
764    #[must_use]
765    fn inverse(&self) -> Self;
766
767    /// Adjoint representation element (for matrix groups: conjugate transpose).
768    ///
769    /// For matrix groups: `conjugate_transpose(g) = g†`
770    /// For unitary groups: `conjugate_transpose(g) = inverse(g)`
771    ///
772    /// # Gauge Theory Application
773    ///
774    /// In gauge transformations, fields transform as:
775    /// ```text
776    /// A' = g A g† + g dg†
777    /// ```
778    /// where g† is the conjugate transpose.
779    ///
780    /// # Examples
781    ///
782    /// ```
783    /// use lie_groups::{LieGroup, SU2};
784    ///
785    /// let g = SU2::rotation_x(0.7);
786    /// let g_dag = g.conjugate_transpose();
787    ///
788    /// // For unitary groups, conjugate transpose equals inverse
789    /// assert!(g_dag.compose(&g).distance_to_identity() < 1e-10);
790    /// ```
791    #[must_use]
792    fn conjugate_transpose(&self) -> Self;
793
794    /// Adjoint representation: `Ad_g`: 𝔤 → 𝔤
795    ///
796    /// The adjoint representation is a group homomorphism `Ad`: G → Aut(𝔤) that maps
797    /// each group element g to a linear automorphism of the Lie algebra.
798    ///
799    /// # Mathematical Definition
800    ///
801    /// For matrix Lie groups:
802    /// ```text
803    /// Ad_g(X) = g X g⁻¹
804    /// ```
805    ///
806    /// For unitary groups where g† = g⁻¹:
807    /// ```text
808    /// Ad_g(X) = g X g†
809    /// ```
810    ///
811    /// # Distinction from `conjugate_transpose()`
812    ///
813    /// - `conjugate_transpose()`: Returns g† (conjugate transpose), a **group element**
814    /// - `adjoint_action()`: Returns `Ad_g(X)` (conjugation), a **Lie algebra element**
815    ///
816    /// | Method | Signature | Output Type | Mathematical Meaning |
817    /// |--------|-----------|-------------|---------------------|
818    /// | `conjugate_transpose()` | `g → g†` | Group element | Hermitian conjugate |
819    /// | `adjoint_action()` | `(g, X) → gXg⁻¹` | Algebra element | Conjugation action |
820    ///
821    /// # Properties
822    ///
823    /// The adjoint representation is a group homomorphism:
824    /// ```text
825    /// Ad_{g₁g₂}(X) = Ad_{g₁}(Ad_{g₂}(X))
826    /// Ad_e(X) = X
827    /// Ad_g⁻¹(X) = Ad_g⁻¹(X)
828    /// ```
829    ///
830    /// It preserves the Lie bracket:
831    /// ```text
832    /// Ad_g([X, Y]) = [Ad_g(X), Ad_g(Y)]
833    /// ```
834    ///
835    /// For abelian groups (e.g., U(1)):
836    /// ```text
837    /// Ad_g(X) = X  (trivial action)
838    /// ```
839    ///
840    /// # Gauge Theory Application
841    ///
842    /// In gauge transformations, the connection A ∈ 𝔤 transforms as:
843    /// ```text
844    /// A' = Ad_g(A) + g dg†
845    ///    = g A g† + (derivative term)
846    /// ```
847    ///
848    /// The adjoint action captures how gauge fields rotate under symmetry transformations.
849    ///
850    /// # Examples
851    ///
852    /// ```ignore
853    /// use lie_groups::{LieGroup, SU2, Su2Algebra};
854    ///
855    /// let g = SU2::rotation_z(std::f64::consts::PI / 2.0);  // 90° rotation
856    /// let X = Su2Algebra::basis_element(0);  // σ_x direction
857    ///
858    /// // Conjugation rotates the algebra element
859    /// let rotated_X = g.adjoint_action(&X);
860    ///
861    /// // Should now point in σ_y direction (rotated by 90°)
862    /// ```
863    ///
864    /// # Complexity
865    ///
866    /// | Group | Time | Space | Notes |
867    /// |-------|------|-------|-------|
868    /// | SU(2) via quaternion | O(1) | O(1) | Quaternion sandwich: qvq* |
869    /// | SU(N) via matrix | O(N³) | O(N²) | Two matrix multiplications |
870    /// | U(1) | O(1) | O(1) | Trivial (abelian) |
871    /// | SO(3) | O(1) | O(1) | 3×3 matrix rotation |
872    ///
873    /// # See Also
874    ///
875    /// - [`exp`](Self::exp) - Exponential map 𝔤 → G
876    /// - [`conjugate_transpose`](Self::conjugate_transpose) - Conjugate transpose g → g†
877    /// - [`LieAlgebra::bracket`] - Lie bracket [·,·]
878    #[must_use]
879    fn adjoint_action(&self, algebra_element: &Self::Algebra) -> Self::Algebra;
880
881    /// Geodesic distance from identity: d(g, e)
882    ///
883    /// Measures the "size" of the group element as the length of the
884    /// shortest geodesic curve from the identity to g on the Lie group
885    /// manifold.
886    ///
887    /// # Mathematical Definition
888    ///
889    /// For matrix Lie groups, typically computed as:
890    /// ```text
891    /// d(g, e) = ||log(g)||_F
892    /// ```
893    /// where ||·||_F is the Frobenius norm and log is the matrix logarithm.
894    ///
895    /// For SU(2), this simplifies to:
896    /// ```text
897    /// d(g, e) = arccos(Re(Tr(g))/2)
898    /// ```
899    /// which is the rotation angle.
900    ///
901    /// # Properties
902    ///
903    /// - Non-negative: `d(g, e) ≥ 0`
904    /// - Identity: `d(e, e) = 0`
905    /// - Symmetric: `d(g, e) = d(g⁻¹, e)`
906    ///
907    /// # Returns
908    ///
909    /// Distance in radians (for rotation groups) or appropriate units.
910    ///
911    /// # Examples
912    ///
913    /// ```
914    /// use lie_groups::{LieGroup, SU2};
915    ///
916    /// let identity = SU2::identity();
917    /// assert!(identity.distance_to_identity().abs() < 1e-10);
918    ///
919    /// let rotation = SU2::rotation_x(1.0);
920    /// assert!((rotation.distance_to_identity() - 1.0).abs() < 1e-10);
921    /// ```
922    #[must_use]
923    fn distance_to_identity(&self) -> f64;
924
925    /// Distance between two group elements: d(g, h)
926    ///
927    /// Computed as the distance from the identity to g⁻¹h:
928    /// ```text
929    /// d(g, h) = d(g⁻¹h, e)
930    /// ```
931    ///
932    /// This is the canonical left-invariant metric on the Lie group.
933    ///
934    /// # Default Implementation
935    ///
936    /// Provided via left-invariance of the metric. Can be overridden
937    /// for efficiency if a direct formula is available.
938    #[must_use]
939    fn distance(&self, other: &Self) -> f64 {
940        self.inverse().compose(other).distance_to_identity()
941    }
942
943    /// Check if this element is approximately the identity.
944    ///
945    /// Useful for numerical algorithms to check convergence or
946    /// validate gauge fixing.
947    ///
948    /// # Parameters
949    ///
950    /// - `tolerance`: Maximum allowed distance from identity
951    ///
952    /// # Examples
953    ///
954    /// ```
955    /// use lie_groups::{LieGroup, SU2};
956    ///
957    /// let almost_identity = SU2::rotation_x(1e-12);
958    /// assert!(almost_identity.is_near_identity(1e-10));
959    /// ```
960    #[must_use]
961    fn is_near_identity(&self, tolerance: f64) -> bool {
962        self.distance_to_identity() < tolerance
963    }
964
965    /// Exponential map: 𝔤 → G
966    ///
967    /// Maps elements from the Lie algebra (tangent space at identity) to the Lie group
968    /// via the matrix exponential or equivalent operation.
969    ///
970    /// # Mathematical Definition
971    ///
972    /// For matrix groups:
973    /// ```text
974    /// exp(X) = I + X + X²/2! + X³/3! + ...
975    /// ```
976    ///
977    /// For SU(2) with X = θ·n̂·(iσ/2) (rotation by angle θ around axis n̂):
978    /// ```text
979    /// exp(X) = cos(θ/2)I + i·sin(θ/2)·(n̂·σ)
980    /// ```
981    ///
982    /// For U(1) with X = iθ:
983    /// ```text
984    /// exp(iθ) = e^(iθ)
985    /// ```
986    ///
987    /// # Properties
988    ///
989    /// - exp(0) = identity
990    /// - exp is a local diffeomorphism near 0
991    /// - exp(t·X) traces a geodesic on the group manifold
992    ///
993    /// # Applications
994    ///
995    /// - **Optimization**: Gradient descent on manifolds via retraction
996    /// - **Integration**: Solving ODEs on Lie groups
997    /// - **Parameterization**: Unconstrained → constrained optimization
998    ///
999    /// # Complexity
1000    ///
1001    /// | Group | Time | Space | Notes |
1002    /// |-------|------|-------|-------|
1003    /// | SU(2) via quaternion | O(1) | O(1) | Rodrigues formula |
1004    /// | SU(N) via Padé | O(N³) | O(N²) | Matrix exp via scaling+squaring |
1005    /// | U(1) | O(1) | O(1) | e^(iθ) |
1006    /// | SO(3) via Rodrigues | O(1) | O(1) | Closed-form rotation |
1007    ///
1008    /// # Examples
1009    ///
1010    /// ```ignore
1011    /// use lie_groups::{LieGroup, SU2, Su2Algebra};
1012    ///
1013    /// let tangent = Su2Algebra::from_components(&[1.0, 0.0, 0.0]);
1014    /// let group_element = SU2::exp(&tangent);
1015    /// // group_element is a rotation by 1 radian around X-axis
1016    /// ```
1017    #[must_use]
1018    fn exp(tangent: &Self::Algebra) -> Self;
1019
1020    /// Logarithm map: G → 𝔤 (inverse of exponential)
1021    ///
1022    /// Maps group elements near the identity back to the Lie algebra.
1023    /// This is the **inverse** of the exponential map in a neighborhood of the identity.
1024    ///
1025    /// # Mathematical Definition
1026    ///
1027    /// For g ∈ G sufficiently close to the identity, returns X ∈ 𝔤 such that:
1028    /// ```text
1029    /// exp(X) = g
1030    /// ```
1031    ///
1032    /// For matrix groups, this is the matrix logarithm:
1033    /// ```text
1034    /// log(g) = Σ_{n=1}^∞ (-1)^{n+1} (g - I)^n / n
1035    /// ```
1036    ///
1037    /// For specific groups:
1038    ///
1039    /// - **U(1)**: log(e^{iθ}) = iθ (principal branch)
1040    /// - **SU(2)**: log(exp(θ n̂·σ)) = θ n̂·σ for θ ∈ [0, π)
1041    /// - **SO(3)**: Rodrigues' formula in reverse
1042    ///
1043    /// # Domain Restrictions
1044    ///
1045    /// Unlike `exp` which is defined on all of 𝔤, `log` is only well-defined
1046    /// on a **neighborhood of the identity**:
1047    ///
1048    /// - **U(1)**: All elements except negative real axis
1049    /// - **SU(2)**: All elements except -I (rotation by π with ambiguous axis)
1050    /// - **SU(3)**: Elements with all eigenvalues in a small neighborhood
1051    ///
1052    /// # Errors
1053    ///
1054    /// Returns `Err(LogError::NotNearIdentity)` if the element is too far from identity.
1055    /// Returns `Err(LogError::Singularity)` if at an exact singularity (e.g., -I for SU(2)).
1056    ///
1057    /// # Gauge Theory Application
1058    ///
1059    /// In lattice gauge theory, the **curvature** (field strength) is computed
1060    /// from the **holonomy** (Wilson loop) via the logarithm:
1061    ///
1062    /// ```text
1063    /// F_□ = log(U_□) / a²
1064    /// ```
1065    ///
1066    /// where:
1067    /// - U_□ ∈ G is the plaquette holonomy (group element)
1068    /// - F_□ ∈ 𝔤 is the curvature (Lie algebra element)
1069    /// - a is the lattice spacing
1070    ///
1071    /// This is the **discrete analog** of the continuum formula:
1072    /// ```text
1073    /// F_{μν} = ∂_μ A_ν - ∂_ν A_μ + [A_μ, A_ν] ∈ 𝔤
1074    /// ```
1075    ///
1076    /// # Numerical Considerations
1077    ///
1078    /// - For elements very close to identity, use Taylor expansion of log
1079    /// - For elements near singularities, logarithm is numerically unstable
1080    /// - In lattice gauge theory, use smaller lattice spacing if log fails
1081    ///
1082    /// # Complexity
1083    ///
1084    /// | Group | Time | Space | Notes |
1085    /// |-------|------|-------|-------|
1086    /// | SU(2) via quaternion | O(1) | O(1) | atan2 + normalize |
1087    /// | SU(N) via eigendecomp | O(N³) | O(N²) | Schur decomposition |
1088    /// | U(1) | O(1) | O(1) | atan2 |
1089    /// | SO(3) via Rodrigues | O(1) | O(1) | Closed-form |
1090    ///
1091    /// # Properties
1092    ///
1093    /// - log(identity) = 0
1094    /// - log(exp(X)) = X for small X (in the cut)
1095    /// - log(g^{-1}) = -log(g) for g near identity
1096    /// - log is a local diffeomorphism near the identity
1097    ///
1098    /// # Examples
1099    ///
1100    /// ```ignore
1101    /// use lie_groups::{LieGroup, SU2, Su2Algebra};
1102    ///
1103    /// // Create rotation around X-axis by 0.5 radians
1104    /// let tangent = Su2Algebra::from_components(&[0.5, 0.0, 0.0]);
1105    /// let g = SU2::exp(&tangent);
1106    ///
1107    /// // Recover the algebra element
1108    /// let recovered = g.log().unwrap();
1109    ///
1110    /// // Should match original (up to numerical precision)
1111    /// assert!(tangent.add(&recovered.scale(-1.0)).norm() < 1e-10);
1112    /// ```
1113    ///
1114    /// # See Also
1115    ///
1116    /// - [`exp`](Self::exp) - Exponential map 𝔤 → G
1117    /// - [`LogError`](crate::LogError) - Error types for logarithm failures
1118    fn log(&self) -> crate::error::LogResult<Self::Algebra>;
1119
1120    /// Trace of the identity element
1121    ///
1122    /// For SU(N), SO(N), and related matrix groups, the identity matrix I has trace N.
1123    /// This is used in gauge theory computations:
1124    /// - Yang-Mills action: `S_p = N - Re Tr(U_p)` (measures deviation from identity)
1125    /// - Field strength: `|F| ∝ √(N - Re Tr(U))`
1126    ///
1127    /// # Examples
1128    ///
1129    /// ```
1130    /// use lie_groups::{SU2, SU3, LieGroup};
1131    ///
1132    /// assert_eq!(SU2::trace_identity(), 2.0);
1133    /// assert_eq!(SU3::trace_identity(), 3.0);
1134    /// ```
1135    #[must_use]
1136    fn trace_identity() -> f64 {
1137        Self::MATRIX_DIM as f64
1138    }
1139
1140    /// Project element back onto the group manifold using Gram-Schmidt orthogonalization.
1141    ///
1142    /// **Critical for numerical stability** (Tao priority): After many `compose()` or `exp()` operations,
1143    /// floating-point errors accumulate and the result may drift from the group manifold.
1144    /// This method projects the element back onto the manifold, restoring group constraints.
1145    ///
1146    /// # Mathematical Background
1147    ///
1148    /// For matrix Lie groups, the group manifold is defined by constraints:
1149    /// - **SU(N)**: U†U = I (unitary) and det(U) = 1 (special)
1150    /// - **SO(N)**: R†R = I (orthogonal) and det(R) = +1 (special)
1151    ///
1152    /// Gram-Schmidt reorthogonalization ensures U†U = I by:
1153    /// 1. Orthogonalizing columns via modified Gram-Schmidt
1154    /// 2. Rescaling det(U) → 1 if needed
1155    ///
1156    /// # When to Use
1157    ///
1158    /// Call after:
1159    /// - Long chains: `g₁.compose(&g₂).compose(&g₃)....compose(&gₙ)` for large n
1160    /// - Many exponentials: `exp(X₁).compose(&exp(X₂))....`
1161    /// - Scaling-and-squaring: After matrix powers in `exp()` algorithm
1162    ///
1163    /// # Performance vs Accuracy Trade-off
1164    ///
1165    /// - **Cost**: O(N³) for N×N matrices (Gram-Schmidt)
1166    /// - **Benefit**: Prevents catastrophic error accumulation
1167    /// - **Guideline**: Reorthogonalize every ~100 operations in long chains
1168    ///
1169    /// # Example
1170    ///
1171    /// ```ignore
1172    /// use lie_groups::{LieGroup, SU2};
1173    ///
1174    /// // Long composition chain
1175    /// let mut g = SU2::identity();
1176    /// for i in 0..1000 {
1177    ///     g = g.compose(&SU2::rotation_x(0.01));
1178    ///     if i % 100 == 0 {
1179    ///         g = g.reorthogonalize();  // Prevent drift every 100 steps
1180    ///     }
1181    /// }
1182    ///
1183    /// // Verify still on manifold
1184    /// assert!(g.conjugate_transpose().compose(&g).distance_to_identity() < 1e-12);
1185    /// ```
1186    ///
1187    /// # Default Implementation
1188    ///
1189    /// The default implementation uses the exponential/logarithm round-trip:
1190    /// ```text
1191    /// reorthogonalize(g) = exp(log(g))
1192    /// ```
1193    ///
1194    /// This works because:
1195    /// 1. `log(g)` extracts the Lie algebra element (always exact algebraically)
1196    /// 2. `exp(X)` produces an element exactly on the manifold
1197    ///
1198    /// Groups may override this with more efficient methods (e.g., direct Gram-Schmidt for SU(N)).
1199    ///
1200    /// # Error Handling
1201    ///
1202    /// - **Never fails** for elements close to the manifold
1203    /// - For elements far from manifold (catastrophic numerical error), returns best approximation
1204    /// - If `log()` fails (element too far from identity), falls back to identity
1205    ///
1206    /// # See Also
1207    ///
1208    /// - [`exp`](Self::exp) - Exponential map (always produces elements on manifold)
1209    /// - [`log`](Self::log) - Logarithm map (inverse of exp near identity)
1210    /// - [Modified Gram-Schmidt Algorithm](https://en.wikipedia.org/wiki/Gram%E2%80%93Schmidt_process)
1211    ///
1212    /// # References
1213    ///
1214    /// - Tao, "Topics in Random Matrix Theory" (2012) - Numerical stability
1215    /// - Higham, "Accuracy and Stability of Numerical Algorithms" (2002) - Matrix orthogonalization
1216    #[must_use]
1217    fn reorthogonalize(&self) -> Self {
1218        // Default: Round-trip through log-exp
1219        // This works because exp() always produces elements on the manifold
1220        match self.log() {
1221            Ok(algebra) => Self::exp(&algebra),
1222            Err(_) => {
1223                // Fallback: If log fails (element too corrupted), return identity
1224                Self::identity()
1225            }
1226        }
1227    }
1228
1229    /// Geodesic interpolation between two group elements.
1230    ///
1231    /// Computes the point at parameter `t` along the geodesic from `self` to `other`:
1232    /// ```text
1233    /// γ(t) = g · exp(t · log(g⁻¹ · h))
1234    /// ```
1235    ///
1236    /// - `t = 0.0` → returns `self`
1237    /// - `t = 1.0` → returns `other`
1238    /// - `t = 0.5` → geodesic midpoint
1239    ///
1240    /// # Errors
1241    ///
1242    /// Returns `None` if `g⁻¹h` is at the cut locus (log fails).
1243    #[must_use]
1244    fn geodesic(&self, other: &Self, t: f64) -> Option<Self> {
1245        let delta = self.inverse().compose(other);
1246        let tangent = delta.log().ok()?;
1247        Some(self.compose(&Self::exp(&tangent.scale(t))))
1248    }
1249}
1250
1251#[cfg(test)]
1252mod tests {
1253    use super::*;
1254
1255    /// Mock Lie algebra for testing (u(1) ≅ ℝ)
1256    #[derive(Clone, Debug, PartialEq)]
1257    struct U1Algebra {
1258        value: f64,
1259    }
1260
1261    impl LieAlgebra for U1Algebra {
1262        const DIM: usize = 1;
1263
1264        fn zero() -> Self {
1265            Self { value: 0.0 }
1266        }
1267
1268        fn add(&self, other: &Self) -> Self {
1269            Self {
1270                value: self.value + other.value,
1271            }
1272        }
1273
1274        fn scale(&self, scalar: f64) -> Self {
1275            Self {
1276                value: self.value * scalar,
1277            }
1278        }
1279
1280        fn norm(&self) -> f64 {
1281            self.value.abs()
1282        }
1283
1284        fn basis_element(i: usize) -> Self {
1285            assert_eq!(i, 0, "U(1) algebra is 1-dimensional");
1286            Self { value: 1.0 }
1287        }
1288
1289        fn from_components(components: &[f64]) -> Self {
1290            assert_eq!(components.len(), 1);
1291            Self {
1292                value: components[0],
1293            }
1294        }
1295
1296        fn to_components(&self) -> Vec<f64> {
1297            vec![self.value]
1298        }
1299
1300        fn bracket(&self, _other: &Self) -> Self {
1301            Self::zero() // u(1) is abelian
1302        }
1303    }
1304
1305    /// Mock Lie group for testing trait laws
1306    #[derive(Clone, Debug, PartialEq)]
1307    struct U1 {
1308        angle: f64, // U(1) = circle group, element = e^{iθ}
1309    }
1310
1311    impl LieGroup for U1 {
1312        const MATRIX_DIM: usize = 1;
1313
1314        type Algebra = U1Algebra;
1315
1316        fn identity() -> Self {
1317            Self { angle: 0.0 }
1318        }
1319
1320        fn compose(&self, other: &Self) -> Self {
1321            Self {
1322                angle: (self.angle + other.angle) % (2.0 * std::f64::consts::PI),
1323            }
1324        }
1325
1326        fn inverse(&self) -> Self {
1327            Self {
1328                angle: (-self.angle).rem_euclid(2.0 * std::f64::consts::PI),
1329            }
1330        }
1331
1332        fn conjugate_transpose(&self) -> Self {
1333            self.inverse() // U(1) is abelian, so conjugate transpose = inverse
1334        }
1335
1336        fn adjoint_action(&self, algebra_element: &Self::Algebra) -> Self::Algebra {
1337            // Trivial adjoint action for abelian group
1338            algebra_element.clone()
1339        }
1340
1341        fn distance_to_identity(&self) -> f64 {
1342            // Shortest arc distance on circle [0, 2π]
1343            let normalized = self.angle.rem_euclid(2.0 * std::f64::consts::PI);
1344            let dist = normalized.min(2.0 * std::f64::consts::PI - normalized);
1345            dist.abs()
1346        }
1347
1348        fn exp(tangent: &Self::Algebra) -> Self {
1349            // exp(iθ) = e^(iθ), represented as angle θ
1350            Self {
1351                angle: tangent.value.rem_euclid(2.0 * std::f64::consts::PI),
1352            }
1353        }
1354
1355        fn log(&self) -> crate::error::LogResult<Self::Algebra> {
1356            // For U(1), log(e^{iθ}) = iθ
1357            Ok(U1Algebra { value: self.angle })
1358        }
1359    }
1360
1361    #[test]
1362    fn test_identity_law() {
1363        let g = U1 { angle: 1.0 };
1364        let e = U1::identity();
1365
1366        let g_times_e = g.compose(&e);
1367        assert!((g_times_e.angle - g.angle).abs() < 1e-10);
1368    }
1369
1370    #[test]
1371    fn test_inverse_law() {
1372        let g = U1 { angle: 2.5 };
1373        let g_inv = g.inverse();
1374        let product = g.compose(&g_inv);
1375
1376        assert!(product.is_near_identity(1e-10));
1377    }
1378
1379    #[test]
1380    fn test_associativity() {
1381        let g1 = U1 { angle: 0.5 };
1382        let g2 = U1 { angle: 1.2 };
1383        let g3 = U1 { angle: 0.8 };
1384
1385        let left = g1.compose(&g2.compose(&g3));
1386        let right = g1.compose(&g2).compose(&g3);
1387
1388        assert!((left.angle - right.angle).abs() < 1e-10);
1389    }
1390
1391    #[test]
1392    fn test_distance_symmetry() {
1393        let g = U1 { angle: 1.5 };
1394        let g_inv = g.inverse();
1395
1396        let d1 = g.distance_to_identity();
1397        let d2 = g_inv.distance_to_identity();
1398
1399        assert!((d1 - d2).abs() < 1e-10);
1400    }
1401
1402    #[test]
1403    fn test_is_near_identity() {
1404        let e = U1::identity();
1405        assert!(e.is_near_identity(1e-10));
1406
1407        let g = U1 { angle: 1e-12 };
1408        assert!(g.is_near_identity(1e-10));
1409
1410        let h = U1 { angle: 0.1 };
1411        assert!(!h.is_near_identity(1e-10));
1412    }
1413
1414    #[test]
1415    fn test_exponential_map() {
1416        // exp(0) = identity
1417        let zero = U1Algebra::zero();
1418        let exp_zero = U1::exp(&zero);
1419        assert!(exp_zero.is_near_identity(1e-10));
1420
1421        // exp(iθ) = e^(iθ)
1422        let tangent = U1Algebra { value: 1.5 };
1423        let group_elem = U1::exp(&tangent);
1424        assert!((group_elem.angle - 1.5).abs() < 1e-10);
1425    }
1426
1427    #[test]
1428    fn test_matrix_dim() {
1429        assert_eq!(U1::MATRIX_DIM, 1);
1430    }
1431
1432    #[test]
1433    fn test_lie_algebra_zero() {
1434        let zero = U1Algebra::zero();
1435        assert_eq!(zero.value, 0.0);
1436    }
1437
1438    #[test]
1439    fn test_lie_algebra_add() {
1440        let v = U1Algebra { value: 1.0 };
1441        let w = U1Algebra { value: 2.0 };
1442        let sum = v.add(&w);
1443        assert_eq!(sum.value, 3.0);
1444    }
1445
1446    #[test]
1447    fn test_lie_algebra_scale() {
1448        let v = U1Algebra { value: 2.0 };
1449        let scaled = v.scale(3.0);
1450        assert_eq!(scaled.value, 6.0);
1451    }
1452
1453    #[test]
1454    fn test_lie_algebra_norm() {
1455        let v = U1Algebra { value: -3.0 };
1456        assert_eq!(v.norm(), 3.0);
1457    }
1458
1459    #[test]
1460    fn test_lie_algebra_basis() {
1461        let basis = U1Algebra::basis_element(0);
1462        assert_eq!(basis.value, 1.0);
1463    }
1464
1465    #[test]
1466    fn test_lie_algebra_from_components() {
1467        let v = U1Algebra::from_components(&[4.5]);
1468        assert_eq!(v.value, 4.5);
1469    }
1470
1471    #[test]
1472    fn test_matrix_dim_available_at_compile_time() {
1473        assert_eq!(U1::MATRIX_DIM, 1);
1474    }
1475
1476    // Note: Associated const equality (e.g., `LieGroup<DIM = 1>`) is not yet stable
1477    // See https://github.com/rust-lang/rust/issues/92827
1478    //
1479    // Once stabilized, we can write functions like:
1480    //
1481    // ```ignore
1482    // fn dimension_one_only<G: LieGroup<DIM = 1>>() -> usize {
1483    //     G::DIM
1484    // }
1485    // ```
1486    //
1487    // This would enable compile-time verification of dimension compatibility.
1488
1489    // ========================================================================
1490    // Reorthogonalization Tests (TAO PRIORITY - Numerical Stability)
1491    // ========================================================================
1492
1493    /// Test that reorthogonalization fixes drift in long composition chains
1494    #[test]
1495    fn test_reorthogonalize_prevents_drift_su2() {
1496        use crate::SU2;
1497
1498        // Create a long composition chain without reorthogonalization
1499        let mut g_no_reorth = SU2::identity();
1500        for _ in 0..1000 {
1501            g_no_reorth = g_no_reorth.compose(&SU2::rotation_x(0.01));
1502        }
1503
1504        // Check unitarity: U†U should be identity
1505        let unitarity_error_no_reorth = g_no_reorth
1506            .conjugate_transpose()
1507            .compose(&g_no_reorth)
1508            .distance_to_identity();
1509
1510        // Same chain WITH reorthogonalization every 100 steps
1511        let mut g_with_reorth = SU2::identity();
1512        for i in 0..1000 {
1513            g_with_reorth = g_with_reorth.compose(&SU2::rotation_x(0.01));
1514            if i % 100 == 0 {
1515                g_with_reorth = g_with_reorth.reorthogonalize();
1516            }
1517        }
1518
1519        let unitarity_error_with_reorth = g_with_reorth
1520            .conjugate_transpose()
1521            .compose(&g_with_reorth)
1522            .distance_to_identity();
1523
1524        // Reorthogonalization should reduce drift significantly
1525        assert!(
1526            unitarity_error_with_reorth < unitarity_error_no_reorth,
1527            "Reorthogonalization should reduce drift: {} vs {}",
1528            unitarity_error_with_reorth,
1529            unitarity_error_no_reorth
1530        );
1531
1532        // With reorthogonalization, should stay very close to manifold
1533        // Tolerance allows for numerical variation across library versions
1534        assert!(
1535            unitarity_error_with_reorth < 1e-7,
1536            "With reorthogonalization, drift should be minimal: {}",
1537            unitarity_error_with_reorth
1538        );
1539    }
1540
1541    /// Test that `reorthogonalize()` is idempotent
1542    #[test]
1543    fn test_reorthogonalize_idempotent() {
1544        use crate::SU2;
1545        use std::f64::consts::PI;
1546
1547        let g = SU2::rotation_y(PI / 3.0);
1548
1549        let g1 = g.reorthogonalize();
1550        let g2 = g1.reorthogonalize();
1551
1552        // Applying twice should give same result (idempotent)
1553        assert!(
1554            g1.distance(&g2) < 1e-14,
1555            "Reorthogonalization should be idempotent"
1556        );
1557    }
1558
1559    /// Test reorthogonalization preserves group element
1560    #[test]
1561    fn test_reorthogonalize_preserves_element() {
1562        use crate::SU2;
1563        use std::f64::consts::PI;
1564
1565        let g = SU2::rotation_z(PI / 4.0);
1566        let g_reorth = g.reorthogonalize();
1567
1568        // Should be very close to original
1569        assert!(
1570            g.distance(&g_reorth) < 1e-12,
1571            "Reorthogonalization should preserve element"
1572        );
1573
1574        // Reorthogonalized version should be exactly on manifold
1575        let unitarity_error = g_reorth
1576            .conjugate_transpose()
1577            .compose(&g_reorth)
1578            .distance_to_identity();
1579        assert!(
1580            unitarity_error < 1e-14,
1581            "Reorthogonalized element should be on manifold: {}",
1582            unitarity_error
1583        );
1584    }
1585
1586    // ========================================================================
1587    // Operator overloading tests
1588    // ========================================================================
1589
1590    use crate::so3::{So3Algebra, SO3};
1591    use crate::su2::{Su2Algebra, SU2};
1592    use crate::su3::Su3Algebra;
1593
1594    #[test]
1595    fn test_su2_algebra_operators() {
1596        let x = Su2Algebra([1.0, 2.0, 3.0]);
1597        let y = Su2Algebra([0.5, 1.0, 1.5]);
1598
1599        // Add
1600        let sum = x + y;
1601        assert_eq!(sum.0, [1.5, 3.0, 4.5]);
1602
1603        // Sub
1604        let diff = x - y;
1605        assert_eq!(diff.0, [0.5, 1.0, 1.5]);
1606
1607        // Neg
1608        let neg = -x;
1609        assert_eq!(neg.0, [-1.0, -2.0, -3.0]);
1610
1611        // Scalar mul (both directions)
1612        let scaled = x * 2.0;
1613        assert_eq!(scaled.0, [2.0, 4.0, 6.0]);
1614        let scaled2 = 2.0 * x;
1615        assert_eq!(scaled2.0, [2.0, 4.0, 6.0]);
1616    }
1617
1618    #[test]
1619    fn test_su2_group_mul_operator() {
1620        let g = SU2::rotation_x(0.3);
1621        let h = SU2::rotation_y(0.7);
1622
1623        // Operator should match compose
1624        let product_op = &g * &h;
1625        let product_compose = g.compose(&h);
1626        assert_eq!(product_op.matrix(), product_compose.matrix());
1627
1628        // MulAssign
1629        let mut g2 = g.clone();
1630        g2 *= &h;
1631        assert_eq!(g2.matrix(), product_compose.matrix());
1632    }
1633
1634    #[test]
1635    fn test_so3_algebra_operators() {
1636        let x = So3Algebra([1.0, 0.0, 0.0]);
1637        let y = So3Algebra([0.0, 1.0, 0.0]);
1638
1639        let sum = x + y;
1640        assert_eq!(sum.0, [1.0, 1.0, 0.0]);
1641
1642        let scaled = 3.0 * x;
1643        assert_eq!(scaled.0, [3.0, 0.0, 0.0]);
1644
1645        let neg = -y;
1646        assert_eq!(neg.0, [0.0, -1.0, 0.0]);
1647    }
1648
1649    #[test]
1650    fn test_so3_group_mul_operator() {
1651        let r1 = SO3::rotation_x(0.5);
1652        let r2 = SO3::rotation_z(0.3);
1653        let product = &r1 * &r2;
1654        let expected = r1.compose(&r2);
1655        assert!(product.distance(&expected) < 1e-14);
1656    }
1657
1658    #[test]
1659    fn test_u1_algebra_operators() {
1660        let x = crate::U1Algebra(1.5);
1661        let y = crate::U1Algebra(0.5);
1662
1663        assert_eq!((x + y).0, 2.0);
1664        assert_eq!((x - y).0, 1.0);
1665        assert_eq!((-x).0, -1.5);
1666        assert_eq!((x * 3.0).0, 4.5);
1667        assert_eq!((3.0 * x).0, 4.5);
1668    }
1669
1670    #[test]
1671    fn test_u1_group_mul_operator() {
1672        let g = crate::U1::from_angle(1.0);
1673        let h = crate::U1::from_angle(2.0);
1674        let product = &g * &h;
1675        let expected = g.compose(&h);
1676        assert!(product.distance(&expected) < 1e-14);
1677    }
1678
1679    #[test]
1680    fn test_rplus_algebra_operators() {
1681        let x = crate::RPlusAlgebra(2.0);
1682        let y = crate::RPlusAlgebra(0.5);
1683
1684        assert_eq!((x + y).0, 2.5);
1685        assert_eq!((x - y).0, 1.5);
1686        assert_eq!((-x).0, -2.0);
1687        assert_eq!((x * 4.0).0, 8.0);
1688        assert_eq!((4.0 * x).0, 8.0);
1689    }
1690
1691    #[test]
1692    fn test_su3_algebra_operators() {
1693        let x = Su3Algebra([1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]);
1694        let y = Su3Algebra([0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]);
1695
1696        let sum = x + y;
1697        assert_eq!(sum.0[0], 1.0);
1698        assert_eq!(sum.0[1], 1.0);
1699
1700        let scaled = 2.0 * x;
1701        assert_eq!(scaled.0[0], 2.0);
1702    }
1703
1704    #[test]
1705    fn test_sun_algebra_operators() {
1706        use crate::sun::SunAlgebra;
1707
1708        let x = SunAlgebra::<2>::basis_element(0);
1709        let y = SunAlgebra::<2>::basis_element(1);
1710
1711        let sum = &x + &y;
1712        assert_eq!(sum.coefficients, vec![1.0, 1.0, 0.0]);
1713
1714        let diff = x.clone() - y.clone();
1715        assert_eq!(diff.coefficients, vec![1.0, -1.0, 0.0]);
1716
1717        let neg = -x.clone();
1718        assert_eq!(neg.coefficients, vec![-1.0, 0.0, 0.0]);
1719
1720        let scaled = x.clone() * 5.0;
1721        assert_eq!(scaled.coefficients, vec![5.0, 0.0, 0.0]);
1722
1723        let scaled2 = 5.0 * x;
1724        assert_eq!(scaled2.coefficients, vec![5.0, 0.0, 0.0]);
1725    }
1726
1727    #[test]
1728    fn test_sun_group_mul_operator() {
1729        use crate::sun::SUN;
1730
1731        let g = SUN::<3>::identity();
1732        let h = SUN::<3>::identity();
1733        let product = &g * &h;
1734        assert!(product.distance_to_identity() < 1e-14);
1735    }
1736
1737    #[test]
1738    fn test_algebra_operator_consistency_with_trait_methods() {
1739        // Verify that operators produce the same results as trait methods
1740        let x = Su2Algebra([0.3, -0.7, 1.2]);
1741        let y = Su2Algebra([0.5, 0.1, -0.4]);
1742
1743        // operator + should equal LieAlgebra::add
1744        let op_sum = x + y;
1745        let trait_sum = LieAlgebra::add(&x, &y);
1746        assert_eq!(op_sum.0, trait_sum.0);
1747
1748        // operator * scalar should equal LieAlgebra::scale
1749        let op_scaled = x * 2.5;
1750        let trait_scaled = x.scale(2.5);
1751        assert_eq!(op_scaled.0, trait_scaled.0);
1752    }
1753
1754    // ========================================================================
1755    // Error Path Tests
1756    // ========================================================================
1757
1758    #[test]
1759    fn test_su2_log_returns_err_for_drifted_matrix() {
1760        use crate::SU2;
1761        use num_complex::Complex64;
1762
1763        // Construct a matrix that violates unitarity (drift simulation)
1764        let mut bad = SU2::identity();
1765        bad.matrix[[0, 0]] = Complex64::new(2.0, 0.0); // |cos(θ/2)| > 1
1766        let result = SU2::log(&bad);
1767        assert!(result.is_err(), "log should fail for non-unitary matrix");
1768    }
1769
1770    #[test]
1771    fn test_su2_log_at_near_2pi_returns_err() {
1772        use crate::SU2;
1773
1774        // θ = 2π is the cut locus where sin(θ/2) = 0
1775        let near_cut = SU2::exp(&crate::Su2Algebra([
1776            0.0,
1777            0.0,
1778            std::f64::consts::PI * 2.0 - 1e-14,
1779        ]));
1780        // This may or may not error depending on precision, but should not panic
1781        let _result = SU2::log(&near_cut);
1782    }
1783
1784    #[test]
1785    fn test_so3_log_at_pi_returns_err() {
1786        use crate::so3::SO3;
1787
1788        // 180° rotation is a singularity for SO(3) log (axis ambiguous)
1789        let r = SO3::rotation_x(std::f64::consts::PI);
1790        let result = SO3::log(&r);
1791        assert!(
1792            result.is_err(),
1793            "SO(3) log at π rotation should return singularity error"
1794        );
1795    }
1796
1797    #[test]
1798    fn test_bch_error_for_invalid_order() {
1799        use crate::bch::{bch_safe, BchError};
1800        use crate::Su2Algebra;
1801
1802        let x = Su2Algebra([0.1, 0.0, 0.0]);
1803        let y = Su2Algebra([0.0, 0.1, 0.0]);
1804        let result = bch_safe::<SU2>(&x, &y, 7);
1805        assert!(
1806            matches!(result, Err(BchError::InvalidOrder(7))),
1807            "BCH should reject invalid order"
1808        );
1809    }
1810
1811    // ========================================================================
1812    // Cross-Group Correspondence Tests
1813    // ========================================================================
1814
1815    #[test]
1816    fn test_su2_so3_double_cover() {
1817        // SU(2) is a double cover of SO(3): the same rotation angle
1818        // should produce the same adjoint action on algebra elements.
1819        //
1820        // For any g ∈ SU(2), the map Ad_g: su(2) → su(2) factors through
1821        // the isomorphism su(2) ≅ so(3), giving a rotation in SO(3).
1822        use crate::so3::{So3Algebra, SO3};
1823        use crate::su2::{Su2Algebra, SU2};
1824
1825        let angles = [0.3, 1.0, std::f64::consts::PI / 2.0, 2.5];
1826
1827        for &angle in &angles {
1828            // Same physical rotation via SU(2) and SO(3)
1829            let g_su2 = SU2::rotation_z(angle);
1830            let r_so3 = SO3::rotation_z(angle);
1831
1832            // Act on basis vectors via adjoint
1833            let x_su2 = Su2Algebra([1.0, 0.0, 0.0]);
1834            let x_so3 = So3Algebra([1.0, 0.0, 0.0]);
1835
1836            let ad_su2 = g_su2.adjoint_action(&x_su2);
1837            let ad_so3 = r_so3.adjoint_action(&x_so3);
1838
1839            // The su(2) → so(3) isomorphism maps e_a ↦ -L_a,
1840            // so Ad_SU2 and Ad_SO3 should give the same rotation
1841            // up to the sign convention. Compare magnitudes.
1842            let su2_rotated = [ad_su2.0[0], ad_su2.0[1], ad_su2.0[2]];
1843            let so3_rotated = [ad_so3.0[0], ad_so3.0[1], ad_so3.0[2]];
1844
1845            // The components should match (both rotate x-axis by angle around z)
1846            // su(2) has -ε structure constants, so Ad picks up a sign flip
1847            // relative to so(3). The rotation matrix is the same.
1848            for i in 0..3 {
1849                assert!(
1850                    (su2_rotated[i].abs() - so3_rotated[i].abs()) < 1e-10
1851                        || (su2_rotated[i] + so3_rotated[i]).abs() < 1e-10
1852                        || (su2_rotated[i] - so3_rotated[i]).abs() < 1e-10,
1853                    "SU(2)/SO(3) adjoint mismatch at angle {}: su2={:?} so3={:?}",
1854                    angle,
1855                    su2_rotated,
1856                    so3_rotated
1857                );
1858            }
1859        }
1860    }
1861
1862    #[test]
1863    fn test_su2_minus_identity_is_same_so3_rotation() {
1864        // -I ∈ SU(2) maps to I ∈ SO(3) (kernel of the double cover)
1865        use crate::su2::{Su2Algebra, SU2};
1866
1867        // exp(2π ê₃) = -I in SU(2)
1868        let g = SU2::exp(&Su2Algebra([0.0, 0.0, std::f64::consts::PI]));
1869        let minus_g = SU2::exp(&Su2Algebra([0.0, 0.0, -std::f64::consts::PI]));
1870
1871        // Both should give the same adjoint action (since -I acts trivially)
1872        let x = Su2Algebra([1.0, 2.0, 3.0]);
1873        let ad_g = g.adjoint_action(&x);
1874        let ad_minus_g = minus_g.adjoint_action(&x);
1875
1876        for i in 0..3 {
1877            assert!(
1878                (ad_g.0[i] - ad_minus_g.0[i]).abs() < 1e-10,
1879                "±I should give same adjoint action"
1880            );
1881        }
1882    }
1883}