lie_groups/su2.rs
1//! SU(2): The Special Unitary Group in 2 Dimensions
2//!
3//! This module implements SU(2) group elements using **2×2 complex matrices**.
4//! We define the Pauli generators directly (no external dependencies needed).
5//!
6//! # SU(2) Implementations
7//!
8//! This library provides multiple implementations of SU(2):
9//!
10//! ## Matrix-Based (this module)
11//!
12//! - **`SU2`**: Specialized 2×2 complex matrix implementation
13//! - Uses `ndarray` for matrix operations
14//! - Has convenient constructors (`rotation_x`, `rotation_y`, `rotation_z`)
15//! - Good for learning, verification, and physics applications
16//!
17//! - **`SU2Generic`** ([`crate::sun::SUN<2>`]): Generic N×N matrix implementation
18//! - Part of the unified SU(N) framework
19//! - Same performance as `SU2` (both use matrix multiplication)
20//! - Useful when writing code generic over SU(N)
21//!
22//! ## Quaternion-Based
23//!
24//! - **[`UnitQuaternion`](crate::UnitQuaternion)**: Quaternion representation
25//! - Uses SU(2) ≅ S³ isomorphism (4 real numbers)
26//! - **~3× faster** for multiplication than matrix representation
27//! - Better numerical stability (no matrix decomposition)
28//! - Recommended for performance-critical rotation computations
29//!
30//! # Performance Comparison
31//!
32//! Benchmarked on typical hardware (multiplication, 1000 ops):
33//! ```text
34//! UnitQuaternion: ~16 µs (4 real multiplies + adds)
35//! SU2 (matrix): ~45 µs (complex 2×2 matrix multiply)
36//! ```
37//!
38//! For applications requiring compatibility with higher SU(N) groups (e.g.,
39//! Wilson loops, parallel transport), the matrix representation is preferred.
40//!
41//! ## Example
42//!
43//! ```rust
44//! use lie_groups::{SU2, LieGroup};
45//!
46//! // SU2 has specialized constructors:
47//! let g1 = SU2::rotation_x(0.5);
48//! let g2 = SU2::rotation_y(0.3);
49//!
50//! // Compose rotations
51//! let g3 = g1.compose(&g2);
52//! ```
53
54use crate::traits::{AntiHermitianByConstruction, LieAlgebra, LieGroup, TracelessByConstruction};
55use ndarray::Array2;
56use num_complex::Complex64;
57use std::fmt;
58use std::ops::{Add, Mul, MulAssign, Neg, Sub};
59
60// ============================================================================
61// Numerical Tolerance Constants
62// ============================================================================
63//
64// These thresholds are chosen based on IEEE 754 f64 precision (ε ≈ 2.2e-16)
65// and practical considerations for accumulated floating-point error.
66//
67// For SU(2) matrices, unitarity ensures |Re(Tr(U))/2| ≤ 1 exactly.
68// Violations indicate matrix corruption (drift from the group manifold).
69//
70// Reference: Higham, "Accuracy and Stability of Numerical Algorithms" (2002)
71
72/// Threshold for detecting unitarity violations.
73/// If |cos(θ/2)| exceeds 1 by more than this, the matrix has drifted from SU(2).
74const UNITARITY_VIOLATION_THRESHOLD: f64 = 1e-10;
75
76/// Threshold for small angle detection in `exp()`.
77/// Below this, return identity to avoid division by near-zero.
78const SMALL_ANGLE_EXP_THRESHOLD: f64 = 1e-10;
79
80/// Threshold for sin(θ/2) in `log()` axis extraction.
81/// Below this, the rotation axis cannot be reliably determined.
82const SIN_HALF_THETA_THRESHOLD: f64 = 1e-10;
83
84/// Lie algebra su(2) ≅ ℝ³
85///
86/// The Lie algebra of SU(2) consists of 2×2 traceless anti-Hermitian matrices.
87/// We represent these using 3 real coordinates.
88///
89/// # Convention
90///
91/// We identify su(2) with ℝ³ via the basis `{e₀, e₁, e₂} = {iσₓ/2, iσᵧ/2, iσᵤ/2}`,
92/// where σᵢ are the Pauli matrices:
93/// ```text
94/// σₓ = [[0, 1], [1, 0]] e₀ = iσₓ/2 = [[0, i/2], [i/2, 0]]
95/// σᵧ = [[0, -i], [i, 0]] e₁ = iσᵧ/2 = [[0, 1/2], [-1/2, 0]]
96/// σᵤ = [[1, 0], [0, -1]] e₂ = iσᵤ/2 = [[i/2, 0], [0, -i/2]]
97/// ```
98///
99/// An element `Su2Algebra::new([a, b, c])` corresponds to the matrix
100/// `(a·iσₓ + b·iσᵧ + c·iσᵤ)/2`, and the parameter `‖(a,b,c)‖` is the
101/// rotation angle in the exponential map.
102///
103/// # Structure Constants
104///
105/// With this basis, the Lie bracket satisfies `[eᵢ, eⱼ] = -εᵢⱼₖ eₖ`, giving
106/// structure constants `fᵢⱼₖ = -εᵢⱼₖ` (negative Levi-Civita symbol).
107/// In ℝ³ coordinates, `[X, Y] = -(X × Y)`.
108///
109/// # Isomorphism with ℝ³
110///
111/// su(2) is isomorphic to ℝ³ as a vector space, and as a Lie algebra
112/// the bracket is the negative cross product. The norm `‖v‖` equals the
113/// rotation angle θ, matching the exponential map
114/// `exp(v) = cos(θ/2)I + i·sin(θ/2)·v̂·σ`.
115///
116/// # Examples
117///
118/// ```
119/// use lie_groups::su2::Su2Algebra;
120/// use lie_groups::traits::LieAlgebra;
121///
122/// // Create algebra element in X direction
123/// let v = Su2Algebra::from_components(&[1.0, 0.0, 0.0]);
124///
125/// // Scale and add
126/// let w = v.scale(2.0);
127/// let sum = v.add(&w);
128/// ```
129#[derive(Clone, Copy, Debug, PartialEq)]
130pub struct Su2Algebra(pub(crate) [f64; 3]);
131
132impl Su2Algebra {
133 /// Create a new su(2) algebra element from components.
134 ///
135 /// The components `[a, b, c]` correspond to the element
136 /// `(a·iσₓ + b·iσᵧ + c·iσᵤ)/2` in the Pauli basis.
137 #[must_use]
138 pub fn new(components: [f64; 3]) -> Self {
139 Self(components)
140 }
141
142 /// Returns the components as a fixed-size array reference.
143 #[must_use]
144 pub fn components(&self) -> &[f64; 3] {
145 &self.0
146 }
147}
148
149impl Add for Su2Algebra {
150 type Output = Self;
151 fn add(self, rhs: Self) -> Self {
152 Self([
153 self.0[0] + rhs.0[0],
154 self.0[1] + rhs.0[1],
155 self.0[2] + rhs.0[2],
156 ])
157 }
158}
159
160impl Add<&Su2Algebra> for Su2Algebra {
161 type Output = Su2Algebra;
162 fn add(self, rhs: &Su2Algebra) -> Su2Algebra {
163 Self([
164 self.0[0] + rhs.0[0],
165 self.0[1] + rhs.0[1],
166 self.0[2] + rhs.0[2],
167 ])
168 }
169}
170
171impl Add<Su2Algebra> for &Su2Algebra {
172 type Output = Su2Algebra;
173 fn add(self, rhs: Su2Algebra) -> Su2Algebra {
174 Su2Algebra([
175 self.0[0] + rhs.0[0],
176 self.0[1] + rhs.0[1],
177 self.0[2] + rhs.0[2],
178 ])
179 }
180}
181
182impl Add<&Su2Algebra> for &Su2Algebra {
183 type Output = Su2Algebra;
184 fn add(self, rhs: &Su2Algebra) -> Su2Algebra {
185 Su2Algebra([
186 self.0[0] + rhs.0[0],
187 self.0[1] + rhs.0[1],
188 self.0[2] + rhs.0[2],
189 ])
190 }
191}
192
193impl Sub for Su2Algebra {
194 type Output = Self;
195 fn sub(self, rhs: Self) -> Self {
196 Self([
197 self.0[0] - rhs.0[0],
198 self.0[1] - rhs.0[1],
199 self.0[2] - rhs.0[2],
200 ])
201 }
202}
203
204impl Neg for Su2Algebra {
205 type Output = Self;
206 fn neg(self) -> Self {
207 Self([-self.0[0], -self.0[1], -self.0[2]])
208 }
209}
210
211impl Mul<f64> for Su2Algebra {
212 type Output = Self;
213 fn mul(self, scalar: f64) -> Self {
214 Self([self.0[0] * scalar, self.0[1] * scalar, self.0[2] * scalar])
215 }
216}
217
218impl Mul<Su2Algebra> for f64 {
219 type Output = Su2Algebra;
220 fn mul(self, rhs: Su2Algebra) -> Su2Algebra {
221 rhs * self
222 }
223}
224
225impl LieAlgebra for Su2Algebra {
226 const DIM: usize = 3;
227
228 #[inline]
229 fn zero() -> Self {
230 Self([0.0, 0.0, 0.0])
231 }
232
233 #[inline]
234 fn add(&self, other: &Self) -> Self {
235 Self([
236 self.0[0] + other.0[0],
237 self.0[1] + other.0[1],
238 self.0[2] + other.0[2],
239 ])
240 }
241
242 #[inline]
243 fn scale(&self, scalar: f64) -> Self {
244 Self([self.0[0] * scalar, self.0[1] * scalar, self.0[2] * scalar])
245 }
246
247 #[inline]
248 fn norm(&self) -> f64 {
249 (self.0[0].powi(2) + self.0[1].powi(2) + self.0[2].powi(2)).sqrt()
250 }
251
252 #[inline]
253 fn basis_element(i: usize) -> Self {
254 assert!(i < 3, "SU(2) algebra is 3-dimensional");
255 let mut v = [0.0; 3];
256 v[i] = 1.0;
257 Self(v)
258 }
259
260 #[inline]
261 fn from_components(components: &[f64]) -> Self {
262 assert_eq!(components.len(), 3, "su(2) has dimension 3");
263 Self([components[0], components[1], components[2]])
264 }
265
266 #[inline]
267 fn to_components(&self) -> Vec<f64> {
268 self.0.to_vec()
269 }
270
271 /// Lie bracket for su(2): [X, Y] = -(X × Y)
272 ///
273 /// # Convention
274 ///
275 /// We represent su(2) as ℝ³ with basis `{eᵢ} = {iσᵢ/2}`. The matrix
276 /// commutator gives:
277 /// ```text
278 /// [iσᵢ/2, iσⱼ/2] = (i²/4)[σᵢ, σⱼ] = -(1/4)(2iεᵢⱼₖσₖ) = -εᵢⱼₖ(iσₖ/2)
279 /// ```
280 ///
281 /// In ℝ³ coordinates, this is the **negative cross product**:
282 /// ```text
283 /// [X, Y] = -(X × Y)
284 /// ```
285 ///
286 /// The negative sign is the unique bracket consistent with the half-angle
287 /// exponential map `exp(v) = cos(‖v‖/2)I + i·sin(‖v‖/2)·v̂·σ`, ensuring
288 /// the BCH formula `exp(X)·exp(Y) = exp(X + Y - ½(X×Y) + ...)` holds.
289 ///
290 /// # Properties
291 ///
292 /// - Structure constants: `fᵢⱼₖ = -εᵢⱼₖ`
293 /// - Antisymmetric: `[X, Y] = -[Y, X]`
294 /// - Jacobi identity: `[X, [Y, Z]] + [Y, [Z, X]] + [Z, [X, Y]] = 0`
295 /// - Killing form: `B(X, Y) = -2(X · Y)`
296 ///
297 /// # Examples
298 ///
299 /// ```
300 /// use lie_groups::Su2Algebra;
301 /// use lie_groups::LieAlgebra;
302 ///
303 /// let e1 = Su2Algebra::basis_element(0); // (1, 0, 0)
304 /// let e2 = Su2Algebra::basis_element(1); // (0, 1, 0)
305 /// let bracket = e1.bracket(&e2); // [e₁, e₂] = -e₃
306 ///
307 /// // Should give -e₃ = (0, 0, -1)
308 /// assert!((bracket.components()[0]).abs() < 1e-10);
309 /// assert!((bracket.components()[1]).abs() < 1e-10);
310 /// assert!((bracket.components()[2] - (-1.0)).abs() < 1e-10);
311 /// ```
312 #[inline]
313 fn bracket(&self, other: &Self) -> Self {
314 // Negative cross product: -(X × Y)
315 // Consistent with exp convention exp(v) = exp_matrix(iv·σ/2)
316 let x = self.0;
317 let y = other.0;
318 Self([
319 -(x[1] * y[2] - x[2] * y[1]),
320 -(x[2] * y[0] - x[0] * y[2]),
321 -(x[0] * y[1] - x[1] * y[0]),
322 ])
323 }
324
325 #[inline]
326 fn inner(&self, other: &Self) -> f64 {
327 self.0[0] * other.0[0] + self.0[1] * other.0[1] + self.0[2] * other.0[2]
328 }
329}
330
331// ============================================================================
332// Casimir Operators for SU(2)
333// ============================================================================
334
335impl crate::Casimir for Su2Algebra {
336 type Representation = crate::representation::Spin;
337
338 #[inline]
339 fn quadratic_casimir_eigenvalue(irrep: &Self::Representation) -> f64 {
340 let j = irrep.value();
341 j * (j + 1.0)
342 }
343
344 #[inline]
345 fn rank() -> usize {
346 1 // SU(2) has rank 1 (dimension of Cartan subalgebra)
347 }
348}
349
350/// SU(2) group element represented as a 2×2 unitary matrix.
351///
352/// SU(2) is the group of 2×2 complex unitary matrices with determinant 1:
353/// ```text
354/// U ∈ SU(2) ⟺ U†U = I and det(U) = 1
355/// ```
356///
357/// # Applications
358/// - Represents rotations and spin transformations
359/// - Acts on spinor fields: ψ → U ψ
360/// - Preserves inner products (unitarity)
361#[derive(Debug, Clone, PartialEq)]
362pub struct SU2 {
363 /// The 2×2 unitary matrix representation
364 pub(crate) matrix: Array2<Complex64>,
365}
366
367impl SU2 {
368 /// Access the underlying 2×2 unitary matrix
369 #[must_use]
370 pub fn matrix(&self) -> &Array2<Complex64> {
371 &self.matrix
372 }
373
374 /// Identity element: I₂
375 #[must_use]
376 pub fn identity() -> Self {
377 Self {
378 matrix: Array2::eye(2),
379 }
380 }
381
382 /// Get Pauli X generator (`σ_x` / 2)
383 ///
384 /// Returns: (1/2) * [[0, 1], [1, 0]]
385 #[must_use]
386 pub fn pauli_x() -> Array2<Complex64> {
387 let mut matrix = Array2::zeros((2, 2));
388 matrix[[0, 1]] = Complex64::new(0.5, 0.0);
389 matrix[[1, 0]] = Complex64::new(0.5, 0.0);
390 matrix
391 }
392
393 /// Get Pauli Y generator (`σ_y` / 2)
394 ///
395 /// Returns: (1/2) * [[0, -i], [i, 0]]
396 #[must_use]
397 pub fn pauli_y() -> Array2<Complex64> {
398 let mut matrix = Array2::zeros((2, 2));
399 matrix[[0, 1]] = Complex64::new(0.0, -0.5);
400 matrix[[1, 0]] = Complex64::new(0.0, 0.5);
401 matrix
402 }
403
404 /// Get Pauli Z generator (`σ_z` / 2)
405 ///
406 /// Returns: (1/2) * [[1, 0], [0, -1]]
407 #[must_use]
408 pub fn pauli_z() -> Array2<Complex64> {
409 let mut matrix = Array2::zeros((2, 2));
410 matrix[[0, 0]] = Complex64::new(0.5, 0.0);
411 matrix[[1, 1]] = Complex64::new(-0.5, 0.0);
412 matrix
413 }
414
415 /// Rotation around X-axis by angle θ
416 ///
417 /// Uses closed form: U = cos(θ/2) I + i sin(θ/2) `σ_x`
418 #[inline]
419 #[must_use]
420 pub fn rotation_x(theta: f64) -> Self {
421 let (s, c) = (theta / 2.0).sin_cos();
422
423 let mut matrix = Array2::zeros((2, 2));
424 matrix[[0, 0]] = Complex64::new(c, 0.0);
425 matrix[[0, 1]] = Complex64::new(0.0, s);
426 matrix[[1, 0]] = Complex64::new(0.0, s);
427 matrix[[1, 1]] = Complex64::new(c, 0.0);
428
429 Self { matrix }
430 }
431
432 /// Rotation around Y-axis by angle θ
433 ///
434 /// Uses closed form: U = cos(θ/2) I + i sin(θ/2) `σ_y`
435 #[inline]
436 #[must_use]
437 pub fn rotation_y(theta: f64) -> Self {
438 let (s, c) = (theta / 2.0).sin_cos();
439
440 let mut matrix = Array2::zeros((2, 2));
441 matrix[[0, 0]] = Complex64::new(c, 0.0);
442 matrix[[0, 1]] = Complex64::new(s, 0.0);
443 matrix[[1, 0]] = Complex64::new(-s, 0.0);
444 matrix[[1, 1]] = Complex64::new(c, 0.0);
445
446 Self { matrix }
447 }
448
449 /// Rotation around Z-axis by angle θ
450 ///
451 /// Uses closed form: U = cos(θ/2) I + i sin(θ/2) `σ_z`
452 #[inline]
453 #[must_use]
454 pub fn rotation_z(theta: f64) -> Self {
455 let (s, c) = (theta / 2.0).sin_cos();
456
457 let mut matrix = Array2::zeros((2, 2));
458 matrix[[0, 0]] = Complex64::new(c, s);
459 matrix[[0, 1]] = Complex64::new(0.0, 0.0);
460 matrix[[1, 0]] = Complex64::new(0.0, 0.0);
461 matrix[[1, 1]] = Complex64::new(c, -s);
462
463 Self { matrix }
464 }
465
466 /// Random SU(2) element uniformly distributed according to Haar measure
467 ///
468 /// Requires the `rand` feature (enabled by default).
469 ///
470 /// Samples uniformly from SU(2) ≅ S³ (the 3-sphere) using the quaternion
471 /// representation and Gaussian sampling method.
472 ///
473 /// # Mathematical Background
474 ///
475 /// SU(2) is diffeomorphic to the 3-sphere S³ ⊂ ℝ⁴:
476 /// ```text
477 /// SU(2) = {(a,b,c,d) ∈ ℝ⁴ : a² + b² + c² + d² = 1}
478 /// ```
479 ///
480 /// represented as matrices:
481 /// ```text
482 /// U = [[a+ib, c+id ],
483 /// [-c+id, a-ib]]
484 /// ```
485 ///
486 /// ## Haar Measure Sampling
487 ///
488 /// To sample uniformly from S³:
489 /// 1. Generate 4 independent standard Gaussian variables (μ=0, σ=1)
490 /// 2. Normalize to unit length
491 ///
492 /// This gives the uniform (Haar) measure on SU(2) due to the rotational
493 /// invariance of the Gaussian distribution.
494 ///
495 /// # References
496 ///
497 /// - Shoemake, K.: "Uniform Random Rotations" (Graphics Gems III, 1992)
498 /// - Diaconis & Saloff-Coste: "Comparison Theorems for Random Walks on Groups" (1993)
499 ///
500 /// # Examples
501 ///
502 /// ```rust
503 /// use lie_groups::SU2;
504 /// use rand::SeedableRng;
505 ///
506 /// let mut rng = rand::rngs::StdRng::seed_from_u64(42);
507 /// let g = SU2::random_haar(&mut rng);
508 ///
509 /// // Verify it's unitary
510 /// assert!(g.verify_unitarity(1e-10));
511 /// ```
512 #[cfg(feature = "rand")]
513 #[must_use]
514 pub fn random_haar<R: rand::Rng>(rng: &mut R) -> Self {
515 use rand_distr::{Distribution, StandardNormal};
516
517 // Numerical stability constant: minimum acceptable norm before re-sampling
518 // Probability of norm < 1e-10 for 4 Gaussians is ~10^-40, but we guard anyway
519 const MIN_NORM: f64 = 1e-10;
520
521 loop {
522 // Generate 4 independent Gaussian(0,1) random variables
523 let a: f64 = StandardNormal.sample(rng);
524 let b: f64 = StandardNormal.sample(rng);
525 let c: f64 = StandardNormal.sample(rng);
526 let d: f64 = StandardNormal.sample(rng);
527
528 // Normalize to unit length (project onto S³)
529 let norm = (a * a + b * b + c * c + d * d).sqrt();
530
531 // Guard against numerical instability: if norm is too small, re-sample
532 // This prevents division by ~0 which would produce Inf/NaN
533 if norm < MIN_NORM {
534 continue;
535 }
536
537 let a = a / norm;
538 let b = b / norm;
539 let c = c / norm;
540 let d = d / norm;
541
542 // Construct SU(2) matrix from quaternion (a,b,c,d)
543 // U = [[a+ib, c+id ],
544 // [-c+id, a-ib]]
545 let mut matrix = Array2::zeros((2, 2));
546 matrix[[0, 0]] = Complex64::new(a, b);
547 matrix[[0, 1]] = Complex64::new(c, d);
548 matrix[[1, 0]] = Complex64::new(-c, d);
549 matrix[[1, 1]] = Complex64::new(a, -b);
550
551 return Self { matrix };
552 }
553 }
554
555 /// Group inverse: U⁻¹ = U† (conjugate transpose for unitary matrices)
556 #[must_use]
557 pub fn inverse(&self) -> Self {
558 let matrix = self.matrix.t().mapv(|z| z.conj());
559 Self { matrix }
560 }
561
562 /// Hermitian conjugate transpose: U†
563 ///
564 /// For unitary matrices U ∈ SU(2), we have U† = U⁻¹
565 /// This is used in gauge transformations: A' = g A g†
566 #[must_use]
567 pub fn conjugate_transpose(&self) -> Self {
568 self.inverse()
569 }
570
571 /// Act on a 2D vector: ψ → U ψ
572 #[inline]
573 #[must_use]
574 pub fn act_on_vector(&self, v: &[Complex64; 2]) -> [Complex64; 2] {
575 let m = &self.matrix;
576 [
577 m[[0, 0]] * v[0] + m[[0, 1]] * v[1],
578 m[[1, 0]] * v[0] + m[[1, 1]] * v[1],
579 ]
580 }
581
582 /// Trace of the matrix: Tr(U)
583 #[must_use]
584 pub fn trace(&self) -> Complex64 {
585 self.matrix[[0, 0]] + self.matrix[[1, 1]]
586 }
587
588 /// Distance from identity (geodesic distance in Lie group manifold)
589 ///
590 /// Computed as: d(U, I) = ||log(U)||_F (Frobenius norm of logarithm)
591 ///
592 /// # Numerical Behavior
593 ///
594 /// For valid SU(2) matrices, `|Re(Tr(U))/2| ≤ 1` exactly. Small violations
595 /// (within `UNITARITY_VIOLATION_THRESHOLD`) are clamped silently. Larger
596 /// violations trigger a debug assertion, indicating matrix corruption.
597 #[must_use]
598 pub fn distance_to_identity(&self) -> f64 {
599 // For SU(2), use the formula: d = arccos(Re(Tr(U))/2)
600 // This is the angle of rotation
601 let trace_val = self.trace();
602 let raw_cos_half = trace_val.re / 2.0;
603
604 // Detect unitarity violations: for valid SU(2), |cos(θ/2)| ≤ 1
605 debug_assert!(
606 raw_cos_half.abs() <= 1.0 + UNITARITY_VIOLATION_THRESHOLD,
607 "SU(2) unitarity violation: |cos(θ/2)| = {:.10} > 1 + ε. \
608 Matrix has drifted from the group manifold.",
609 raw_cos_half.abs()
610 );
611
612 // Safe clamp for tiny numerical overshoot
613 let cos_half_angle = raw_cos_half.clamp(-1.0, 1.0);
614 2.0 * cos_half_angle.acos()
615 }
616
617 /// Extract rotation angle θ and axis n̂ from SU(2) element.
618 ///
619 /// For U = cos(θ/2)·I + i·sin(θ/2)·(n̂·σ), returns (θ, n̂).
620 ///
621 /// # Returns
622 ///
623 /// Tuple of:
624 /// - `angle`: Rotation angle θ ∈ [0, 2π]
625 /// - `axis`: Unit 3-vector n̂ = `[n_x, n_y, n_z]` (or `[0,0,1]` if θ ≈ 0)
626 ///
627 /// # Usage in Topological Charge
628 ///
629 /// For computing topological charge, the product `F_μν · F_ρσ` involves
630 /// the dot product of the orientation axes: `(n_μν · n_ρσ)`.
631 #[must_use]
632 pub fn angle_and_axis(&self) -> (f64, [f64; 3]) {
633 // U = cos(θ/2)·I + i·sin(θ/2)·(n̂·σ)
634 // matrix[[0,0]] = cos(θ/2) + i·sin(θ/2)·n_z
635 // matrix[[0,1]] = sin(θ/2)·n_y + i·sin(θ/2)·n_x
636 // matrix[[1,0]] = -sin(θ/2)·n_y + i·sin(θ/2)·n_x
637 // matrix[[1,1]] = cos(θ/2) - i·sin(θ/2)·n_z
638
639 let cos_half = self.matrix[[0, 0]].re.clamp(-1.0, 1.0);
640 let angle = 2.0 * cos_half.acos();
641
642 // For small angles, axis is ill-defined; return default
643 if angle < 1e-10 {
644 return (angle, [0.0, 0.0, 1.0]);
645 }
646
647 let sin_half = (angle / 2.0).sin();
648 if sin_half.abs() < 1e-10 {
649 // Near identity or near -I
650 return (angle, [0.0, 0.0, 1.0]);
651 }
652
653 // Extract axis components
654 let n_z = self.matrix[[0, 0]].im / sin_half;
655 let n_x = self.matrix[[0, 1]].im / sin_half;
656 let n_y = self.matrix[[0, 1]].re / sin_half;
657
658 // Normalize (should already be unit, but ensure numerical stability)
659 let norm = (n_x * n_x + n_y * n_y + n_z * n_z).sqrt();
660 if norm < 1e-10 {
661 return (angle, [0.0, 0.0, 1.0]);
662 }
663
664 (angle, [n_x / norm, n_y / norm, n_z / norm])
665 }
666
667 /// Verify this is approximately in SU(2)
668 ///
669 /// Checks both:
670 /// - Unitarity: U†U ≈ I
671 /// - Special: |det(U) - 1| < tolerance
672 #[must_use]
673 pub fn verify_unitarity(&self, tolerance: f64) -> bool {
674 // Check U†U = I
675 let u_dagger = self.matrix.t().mapv(|z: Complex64| z.conj());
676 let product = u_dagger.dot(&self.matrix);
677 let identity = Array2::eye(2);
678
679 let diff_norm: f64 = product
680 .iter()
681 .zip(identity.iter())
682 .map(|(a, b): (&Complex64, &Complex64)| (a - b).norm_sqr())
683 .sum::<f64>()
684 .sqrt();
685
686 if diff_norm >= tolerance {
687 return false;
688 }
689
690 // Check det(U) = 1 (not just unitary, but special unitary)
691 let det =
692 self.matrix[[0, 0]] * self.matrix[[1, 1]] - self.matrix[[0, 1]] * self.matrix[[1, 0]];
693 (det - Complex64::new(1.0, 0.0)).norm() < tolerance
694 }
695
696 /// Convert to 2×2 array format
697 #[must_use]
698 pub fn to_matrix(&self) -> [[Complex64; 2]; 2] {
699 [
700 [self.matrix[[0, 0]], self.matrix[[0, 1]]],
701 [self.matrix[[1, 0]], self.matrix[[1, 1]]],
702 ]
703 }
704
705 /// Create from 2×2 array format
706 #[must_use]
707 pub fn from_matrix(arr: [[Complex64; 2]; 2]) -> Self {
708 let mut matrix = Array2::zeros((2, 2));
709 matrix[[0, 0]] = arr[0][0];
710 matrix[[0, 1]] = arr[0][1];
711 matrix[[1, 0]] = arr[1][0];
712 matrix[[1, 1]] = arr[1][1];
713 Self { matrix }
714 }
715
716 // ========================================================================
717 // Numerically Stable Logarithm
718 // ========================================================================
719
720 /// Compute logarithm using numerically stable atan2 formulation.
721 ///
722 /// This is more stable than the standard `log()` method, especially for
723 /// rotation angles approaching π. Uses `atan2(sin(θ/2), cos(θ/2))` instead
724 /// of `acos(cos(θ/2))` for improved numerical conditioning.
725 ///
726 /// # Stability Improvements
727 ///
728 /// | Angle θ | acos stability | atan2 stability |
729 /// |---------|----------------|-----------------|
730 /// | θ ≈ 0 | Poor (derivative ∞) | Good |
731 /// | θ ≈ π/2 | Good | Good |
732 /// | θ ≈ π | Good | Good |
733 ///
734 /// The atan2 formulation avoids the derivative singularity at θ = 0
735 /// where d(acos)/dx → ∞.
736 ///
737 /// # Returns
738 ///
739 /// Same as `log()`, but with improved numerical stability.
740 pub fn log_stable(&self) -> crate::error::LogResult<Su2Algebra> {
741 use crate::error::LogError;
742
743 // For SU(2): U = cos(θ/2)I + i·sin(θ/2)·(n̂·σ)
744 //
745 // Extract sin(θ/2)·n̂ directly from matrix elements:
746 // U[[0,0]] = cos(θ/2) + i·nz·sin(θ/2) → sin_nz = Im(U[[0,0]])
747 // U[[0,1]] = (ny + i·nx)·sin(θ/2) → sin_nx = Im(U[[0,1]])
748 // sin_ny = Re(U[[0,1]])
749
750 let sin_nx = self.matrix[[0, 1]].im;
751 let sin_ny = self.matrix[[0, 1]].re;
752 let sin_nz = self.matrix[[0, 0]].im;
753
754 // Compute sin(θ/2) = ||sin(θ/2)·n̂||
755 let sin_half_theta = (sin_nx * sin_nx + sin_ny * sin_ny + sin_nz * sin_nz).sqrt();
756
757 // Compute cos(θ/2) = Re(Tr(U))/2
758 let raw_cos_half = self.matrix[[0, 0]].re; // = cos(θ/2) for proper SU(2)
759
760 // Detect unitarity violations
761 if (raw_cos_half.abs() - 1.0).max(0.0) > UNITARITY_VIOLATION_THRESHOLD
762 && sin_half_theta < UNITARITY_VIOLATION_THRESHOLD
763 {
764 // Check if cos² + sin² ≈ 1
765 let norm_check = raw_cos_half * raw_cos_half + sin_half_theta * sin_half_theta;
766 if (norm_check - 1.0).abs() > 1e-6 {
767 return Err(LogError::NumericalInstability {
768 reason: format!(
769 "SU(2) unitarity violation: cos²(θ/2) + sin²(θ/2) = {:.10} ≠ 1",
770 norm_check
771 ),
772 });
773 }
774 }
775
776 // Use atan2 for stable angle extraction
777 // atan2(sin, cos) is well-conditioned everywhere except at origin
778 let half_theta = sin_half_theta.atan2(raw_cos_half);
779 let theta = 2.0 * half_theta;
780
781 // Handle identity case
782 if sin_half_theta < 1e-15 {
783 return Ok(Su2Algebra::zero());
784 }
785
786 // Note: θ = π is NOT a singularity for SU(2). At θ = π,
787 // sin(θ/2) = 1 and the rotation axis is perfectly extractable.
788 // The true singularity is at θ = 2π (U = -I, sin(θ/2) = 0),
789 // which is caught by the sin_half_theta check above.
790
791 // Extract normalized axis and scale by angle
792 let scale = theta / sin_half_theta;
793 Ok(Su2Algebra([scale * sin_nx, scale * sin_ny, scale * sin_nz]))
794 }
795
796 /// Compute logarithm with conditioning information.
797 ///
798 /// Returns both the logarithm and a [`LogCondition`](crate::LogCondition)
799 /// structure that provides information about the numerical reliability
800 /// of the result.
801 ///
802 /// Unlike `log()`, this method also provides condition number information
803 /// so callers can assess the numerical reliability of the result near
804 /// the cut locus (θ → 2π, U → -I).
805 ///
806 /// # Example
807 ///
808 /// ```rust
809 /// use lie_groups::SU2;
810 ///
811 /// let g = SU2::rotation_x(2.9); // Close to π
812 /// let (log_g, cond) = g.log_with_condition().unwrap();
813 ///
814 /// if cond.is_well_conditioned() {
815 /// println!("Reliable result: {:?}", log_g);
816 /// } else {
817 /// println!("Warning: condition number = {:.1}", cond.condition_number());
818 /// }
819 /// ```
820 ///
821 /// # Cut Locus Behavior
822 ///
823 /// The SU(2) cut locus is at θ = 2π (U = -I), where sin(θ/2) = 0
824 /// and the axis is undefined. As θ → 2π, the condition number
825 /// diverges as 1/sin(θ/2). The method still returns a result — it's
826 /// up to the caller to decide whether to use it based on the
827 /// conditioning information.
828 pub fn log_with_condition(&self) -> crate::error::ConditionedLogResult<Su2Algebra> {
829 use crate::error::{LogCondition, LogError};
830
831 // Extract sin(θ/2)·n̂ from matrix elements
832 let sin_nx = self.matrix[[0, 1]].im;
833 let sin_ny = self.matrix[[0, 1]].re;
834 let sin_nz = self.matrix[[0, 0]].im;
835 let sin_half_theta = (sin_nx * sin_nx + sin_ny * sin_ny + sin_nz * sin_nz).sqrt();
836
837 let raw_cos_half = self.matrix[[0, 0]].re;
838
839 // Unitarity check
840 let norm_check = raw_cos_half * raw_cos_half + sin_half_theta * sin_half_theta;
841 if (norm_check - 1.0).abs() > 1e-6 {
842 return Err(LogError::NumericalInstability {
843 reason: format!(
844 "SU(2) unitarity violation: cos²(θ/2) + sin²(θ/2) = {:.10} ≠ 1",
845 norm_check
846 ),
847 });
848 }
849
850 // Stable angle extraction using atan2
851 let half_theta = sin_half_theta.atan2(raw_cos_half);
852 let theta = 2.0 * half_theta;
853
854 // Compute conditioning
855 let condition = LogCondition::from_angle(theta);
856
857 // Handle identity case
858 if sin_half_theta < 1e-15 {
859 return Ok((Su2Algebra::zero(), condition));
860 }
861
862 // For θ < 2π, the axis is extractable. The result is well-conditioned
863 // for most of the range; only near θ = 2π (U = -I) does sin(θ/2) → 0.
864 let scale = theta / sin_half_theta;
865 let result = Su2Algebra([scale * sin_nx, scale * sin_ny, scale * sin_nz]);
866
867 Ok((result, condition))
868 }
869}
870
871/// Group multiplication: U₁ * U₂
872impl Mul<&SU2> for &SU2 {
873 type Output = SU2;
874
875 fn mul(self, rhs: &SU2) -> SU2 {
876 SU2 {
877 matrix: self.matrix.dot(&rhs.matrix),
878 }
879 }
880}
881
882impl Mul<&SU2> for SU2 {
883 type Output = SU2;
884
885 fn mul(self, rhs: &SU2) -> SU2 {
886 &self * rhs
887 }
888}
889
890impl MulAssign<&SU2> for SU2 {
891 fn mul_assign(&mut self, rhs: &SU2) {
892 self.matrix = self.matrix.dot(&rhs.matrix);
893 }
894}
895
896impl approx::AbsDiffEq for Su2Algebra {
897 type Epsilon = f64;
898
899 fn default_epsilon() -> Self::Epsilon {
900 1e-10
901 }
902
903 fn abs_diff_eq(&self, other: &Self, epsilon: Self::Epsilon) -> bool {
904 self.0
905 .iter()
906 .zip(other.0.iter())
907 .all(|(a, b)| (a - b).abs() < epsilon)
908 }
909}
910
911impl approx::RelativeEq for Su2Algebra {
912 fn default_max_relative() -> Self::Epsilon {
913 1e-10
914 }
915
916 fn relative_eq(
917 &self,
918 other: &Self,
919 epsilon: Self::Epsilon,
920 max_relative: Self::Epsilon,
921 ) -> bool {
922 self.0
923 .iter()
924 .zip(other.0.iter())
925 .all(|(a, b)| approx::RelativeEq::relative_eq(a, b, epsilon, max_relative))
926 }
927}
928
929impl fmt::Display for Su2Algebra {
930 fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
931 write!(
932 f,
933 "su(2)[{:.4}, {:.4}, {:.4}]",
934 self.0[0], self.0[1], self.0[2]
935 )
936 }
937}
938
939impl fmt::Display for SU2 {
940 fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
941 // Show rotation angle and distance to identity
942 let dist = self.distance_to_identity();
943 write!(f, "SU(2)(θ={:.4})", dist)
944 }
945}
946
947/// Implementation of the `LieGroup` trait for SU(2).
948///
949/// This provides the abstract group interface, making SU(2) usable in
950/// generic gauge theory algorithms.
951///
952/// # Mathematical Verification
953///
954/// The implementation satisfies all group axioms:
955/// - Identity: Verified in tests (`test_identity`)
956/// - Inverse: Verified in tests (`test_inverse`)
957/// - Associativity: Follows from matrix multiplication
958/// - Closure: Matrix multiplication of unitaries is unitary
959impl LieGroup for SU2 {
960 const MATRIX_DIM: usize = 2;
961
962 type Algebra = Su2Algebra;
963
964 fn identity() -> Self {
965 Self::identity() // Delegate to inherent method
966 }
967
968 fn compose(&self, other: &Self) -> Self {
969 // Group composition is matrix multiplication
970 self * other
971 }
972
973 fn inverse(&self) -> Self {
974 Self::inverse(self) // Delegate to inherent method
975 }
976
977 fn conjugate_transpose(&self) -> Self {
978 Self::conjugate_transpose(self) // Delegate to inherent method
979 }
980
981 fn distance_to_identity(&self) -> f64 {
982 Self::distance_to_identity(self) // Delegate to inherent method
983 }
984
985 fn exp(tangent: &Su2Algebra) -> Self {
986 // Su2Algebra([a, b, c]) corresponds to X = i(a·σₓ + b·σᵧ + c·σᵤ)/2
987 // with rotation angle θ = ||(a,b,c)||.
988 let angle = tangent.norm();
989
990 if angle < SMALL_ANGLE_EXP_THRESHOLD {
991 // For small angles, return identity (avoid division by zero)
992 return Self::identity();
993 }
994
995 // Normalize to get rotation axis
996 let axis = tangent.scale(1.0 / angle);
997
998 // Rodrigues formula: exp(i(θ/2)·n̂·σ) = cos(θ/2)I + i·sin(θ/2)·n̂·σ
999 let (s, c) = (angle / 2.0).sin_cos();
1000
1001 let mut matrix = Array2::zeros((2, 2));
1002 matrix[[0, 0]] = Complex64::new(c, s * axis.0[2]);
1003 matrix[[0, 1]] = Complex64::new(s * axis.0[1], s * axis.0[0]);
1004 matrix[[1, 0]] = Complex64::new(-s * axis.0[1], s * axis.0[0]);
1005 matrix[[1, 1]] = Complex64::new(c, -s * axis.0[2]);
1006
1007 Self { matrix }
1008 }
1009
1010 fn log(&self) -> crate::error::LogResult<Su2Algebra> {
1011 use crate::error::LogError;
1012
1013 // For SU(2), we use the formula:
1014 // U = cos(θ/2)I + i·sin(θ/2)·(n̂·σ)
1015 //
1016 // Where θ is the rotation angle and n̂ is the rotation axis.
1017 //
1018 // The logarithm is:
1019 // log(U) = θ·n̂ ∈ su(2) ≅ ℝ³
1020 //
1021 // Step 1: Extract angle from trace
1022 // Tr(U) = 2·cos(θ/2)
1023 let tr = self.trace();
1024 let raw_cos_half = tr.re / 2.0;
1025
1026 // Detect unitarity violations: for valid SU(2), |cos(θ/2)| ≤ 1
1027 if raw_cos_half.abs() > 1.0 + UNITARITY_VIOLATION_THRESHOLD {
1028 return Err(LogError::NumericalInstability {
1029 reason: format!(
1030 "SU(2) unitarity violation: |cos(θ/2)| = {:.10} > 1 + ε. \
1031 Matrix has drifted from the group manifold.",
1032 raw_cos_half.abs()
1033 ),
1034 });
1035 }
1036
1037 // Safe clamp for tiny numerical overshoot
1038 let cos_half_theta = raw_cos_half.clamp(-1.0, 1.0);
1039 let half_theta = cos_half_theta.acos();
1040 let theta = 2.0 * half_theta;
1041
1042 // Step 2: Handle special cases
1043 const SMALL_ANGLE_THRESHOLD: f64 = 1e-10;
1044
1045 if theta.abs() < SMALL_ANGLE_THRESHOLD {
1046 // Near identity: log(I) = 0
1047 return Ok(Su2Algebra::zero());
1048 }
1049
1050 // Note: θ = π is NOT a singularity for SU(2). At θ = π,
1051 // sin(θ/2) = 1 and the axis is perfectly extractable.
1052 // The true singularity is at θ = 2π where U = -I and
1053 // sin(θ/2) = 0, which is caught by the sin_half_theta
1054 // check below.
1055
1056 // Step 3: Extract rotation axis from matrix elements
1057 // For U = cos(θ/2)I + i·sin(θ/2)·(n̂·σ), we have:
1058 //
1059 // U = [[cos(θ/2) + i·nz·sin(θ/2), (ny + i·nx)·sin(θ/2) ],
1060 // [(-ny + i·nx)·sin(θ/2), cos(θ/2) - i·nz·sin(θ/2)]]
1061 //
1062 // Extracting:
1063 // U[[0,0]] = cos(θ/2) + i·nz·sin(θ/2) → nz = Im(U[[0,0]]) / sin(θ/2)
1064 // U[[0,1]] = (ny + i·nx)·sin(θ/2) → nx = Im(U[[0,1]]) / sin(θ/2)
1065 // ny = Re(U[[0,1]]) / sin(θ/2)
1066
1067 let sin_half_theta = (half_theta).sin();
1068
1069 if sin_half_theta.abs() < SIN_HALF_THETA_THRESHOLD {
1070 // This shouldn't happen given our checks above, but guard against it
1071 return Err(LogError::NumericalInstability {
1072 reason: "sin(θ/2) too small for reliable axis extraction".to_string(),
1073 });
1074 }
1075
1076 let nx = self.matrix[[0, 1]].im / sin_half_theta;
1077 let ny = self.matrix[[0, 1]].re / sin_half_theta;
1078 let nz = self.matrix[[0, 0]].im / sin_half_theta;
1079
1080 // The logarithm is log(U) = θ·(nx, ny, nz) ∈ su(2)
1081 Ok(Su2Algebra([theta * nx, theta * ny, theta * nz]))
1082 }
1083
1084 fn adjoint_action(&self, algebra_element: &Su2Algebra) -> Su2Algebra {
1085 // Ad_g(X) = g X g† for matrix groups
1086 //
1087 // We construct the matrix M = i(a·σₓ + b·σᵧ + c·σᵤ) = 2X, where
1088 // X = i(a·σₓ + b·σᵧ + c·σᵤ)/2 is the actual algebra element in the
1089 // {iσ/2} basis. The factor of 2 cancels: we compute gMg† and extract
1090 // coefficients using the same convention, recovering (a', b', c').
1091
1092 let [a, b, c] = algebra_element.0;
1093 let i = Complex64::new(0.0, 1.0);
1094
1095 // M = i(a·σₓ + b·σᵧ + c·σᵤ) = [[c·i, b + a·i], [-b + a·i, -c·i]]
1096 let mut x_matrix = Array2::zeros((2, 2));
1097 x_matrix[[0, 0]] = i * c; // c·i
1098 x_matrix[[0, 1]] = Complex64::new(b, a); // b + a·i
1099 x_matrix[[1, 0]] = Complex64::new(-b, a); // -b + a·i
1100 x_matrix[[1, 1]] = -i * c; // -c·i
1101
1102 // Compute g X g† where g† is conjugate transpose
1103 let g_x = self.matrix.dot(&x_matrix);
1104 let g_adjoint_matrix = self.matrix.t().mapv(|z| z.conj()); // Conjugate transpose
1105 let result = g_x.dot(&g_adjoint_matrix);
1106
1107 // Extract coefficients: result = [[c'·i, b'+a'·i], [-b'+a'·i, -c'·i]]
1108
1109 let a_prime = result[[0, 1]].im; // Imaginary part of result[[0,1]]
1110 let b_prime = result[[0, 1]].re; // Real part of result[[0,1]]
1111 let c_prime = result[[0, 0]].im; // Imaginary part of result[[0,0]]
1112
1113 Su2Algebra([a_prime, b_prime, c_prime])
1114 }
1115}
1116
1117// ============================================================================
1118// Mathematical Property Implementations for SU(2)
1119// ============================================================================
1120
1121use crate::traits::{Compact, SemiSimple, Simple};
1122
1123/// SU(2) is compact.
1124///
1125/// All elements are bounded: ||U|| = 1 for all U ∈ SU(2).
1126/// The group is diffeomorphic to the 3-sphere S³.
1127impl Compact for SU2 {}
1128
1129/// SU(2) is simple.
1130///
1131/// It has no non-trivial normal subgroups. This is a fundamental
1132/// result in Lie theory - SU(2) is one of the classical simple groups.
1133impl Simple for SU2 {}
1134
1135/// SU(2) is semi-simple.
1136///
1137/// As a simple group, it is automatically semi-simple.
1138/// (Simple ⊂ Semi-simple in the classification hierarchy)
1139impl SemiSimple for SU2 {}
1140
1141// ============================================================================
1142// Algebra Marker Traits
1143// ============================================================================
1144
1145/// su(2) algebra elements are traceless by construction.
1146///
1147/// The representation `Su2Algebra::new([f64; 3])` stores coefficients in the
1148/// Pauli basis {iσ₁, iσ₂, iσ₃}. Since each Pauli matrix is traceless,
1149/// any linear combination is also traceless.
1150///
1151/// # Lean Connection
1152///
1153/// Combined with `det_exp_eq_exp_trace`: det(exp(X)) = exp(tr(X)) = exp(0) = 1.
1154/// Therefore `SU2::exp` always produces elements with determinant 1.
1155impl TracelessByConstruction for Su2Algebra {}
1156
1157/// su(2) algebra elements are anti-Hermitian by construction.
1158///
1159/// The representation uses {iσ₁, iσ₂, iσ₃} where σᵢ are Hermitian.
1160/// Since (iσ)† = -iσ† = -iσ, each basis element is anti-Hermitian,
1161/// and any real linear combination is also anti-Hermitian.
1162///
1163/// # Lean Connection
1164///
1165/// Combined with `exp_antiHermitian_unitary`: exp(X)† · exp(X) = I.
1166/// Therefore `SU2::exp` always produces unitary elements.
1167impl AntiHermitianByConstruction for Su2Algebra {}
1168
1169#[cfg(test)]
1170mod tests {
1171 use super::*;
1172 use approx::assert_relative_eq;
1173
1174 #[test]
1175 fn test_identity() {
1176 let id = SU2::identity();
1177 assert!(id.verify_unitarity(1e-10));
1178 assert_relative_eq!(id.distance_to_identity(), 0.0, epsilon = 1e-10);
1179 }
1180
1181 #[test]
1182 fn test_rotation_unitarity() {
1183 let u = SU2::rotation_x(0.5);
1184 assert!(u.verify_unitarity(1e-10));
1185 }
1186
1187 #[test]
1188 fn test_inverse() {
1189 let u = SU2::rotation_y(1.2);
1190 let u_inv = u.inverse();
1191 let product = &u * &u_inv;
1192
1193 assert_relative_eq!(product.distance_to_identity(), 0.0, epsilon = 1e-8);
1194 }
1195
1196 #[test]
1197 fn test_group_multiplication() {
1198 let u1 = SU2::rotation_x(0.3);
1199 let u2 = SU2::rotation_y(0.7);
1200 let product = &u1 * &u2;
1201
1202 assert!(product.verify_unitarity(1e-10));
1203 }
1204
1205 #[test]
1206 fn test_action_on_vector() {
1207 let id = SU2::identity();
1208 let v = [Complex64::new(1.0, 0.0), Complex64::new(0.0, 1.0)];
1209 let result = id.act_on_vector(&v);
1210
1211 assert_relative_eq!(result[0].re, v[0].re, epsilon = 1e-10);
1212 assert_relative_eq!(result[1].im, v[1].im, epsilon = 1e-10);
1213 }
1214
1215 #[test]
1216 fn test_rotation_preserves_norm() {
1217 let u = SU2::rotation_z(2.5);
1218 let v = [Complex64::new(3.0, 4.0), Complex64::new(1.0, 2.0)];
1219
1220 let norm_before = v[0].norm_sqr() + v[1].norm_sqr();
1221 let rotated = u.act_on_vector(&v);
1222 let norm_after = rotated[0].norm_sqr() + rotated[1].norm_sqr();
1223
1224 assert_relative_eq!(norm_before, norm_after, epsilon = 1e-10);
1225 }
1226
1227 // ========================================================================
1228 // LieGroup Trait Tests
1229 // ========================================================================
1230
1231 #[test]
1232 fn test_lie_group_identity() {
1233 use crate::traits::LieGroup;
1234
1235 let g = SU2::rotation_x(1.5);
1236 let e = SU2::identity();
1237
1238 // Right identity: g * e = g
1239 let g_times_e = g.compose(&e);
1240 assert_relative_eq!(g_times_e.distance(&g), 0.0, epsilon = 1e-10);
1241
1242 // Left identity: e * g = g
1243 let e_times_g = e.compose(&g);
1244 assert_relative_eq!(e_times_g.distance(&g), 0.0, epsilon = 1e-10);
1245 }
1246
1247 #[test]
1248 fn test_lie_group_inverse() {
1249 use crate::traits::LieGroup;
1250
1251 let g = SU2::rotation_y(2.3);
1252 let g_inv = g.inverse();
1253
1254 // g * g⁻¹ = e
1255 // Use 1e-7 tolerance to account for accumulated numerical error
1256 // in matrix operations (2×2 matrix multiplication + arccos)
1257 let right_product = g.compose(&g_inv);
1258 assert!(
1259 right_product.is_near_identity(1e-7),
1260 "Right inverse failed: distance = {}",
1261 right_product.distance_to_identity()
1262 );
1263
1264 // g⁻¹ * g = e
1265 let left_product = g_inv.compose(&g);
1266 assert!(
1267 left_product.is_near_identity(1e-7),
1268 "Left inverse failed: distance = {}",
1269 left_product.distance_to_identity()
1270 );
1271 }
1272
1273 #[test]
1274 fn test_lie_group_associativity() {
1275 use crate::traits::LieGroup;
1276
1277 let g1 = SU2::rotation_x(0.5);
1278 let g2 = SU2::rotation_y(0.8);
1279 let g3 = SU2::rotation_z(1.2);
1280
1281 // (g1 * g2) * g3
1282 let left_assoc = g1.compose(&g2).compose(&g3);
1283
1284 // g1 * (g2 * g3)
1285 let right_assoc = g1.compose(&g2.compose(&g3));
1286
1287 assert_relative_eq!(left_assoc.distance(&right_assoc), 0.0, epsilon = 1e-10);
1288 }
1289
1290 #[test]
1291 fn test_lie_group_distance_symmetry() {
1292 let g = SU2::rotation_x(1.8);
1293 let g_inv = g.inverse();
1294
1295 // d(g, e) = d(g⁻¹, e) by symmetry
1296 let d1 = g.distance_to_identity();
1297 let d2 = g_inv.distance_to_identity();
1298
1299 assert_relative_eq!(d1, d2, epsilon = 1e-10);
1300 }
1301
1302 #[test]
1303 fn test_lie_group_is_near_identity() {
1304 use crate::traits::LieGroup;
1305
1306 let e = SU2::identity();
1307 assert!(
1308 e.is_near_identity(1e-10),
1309 "Identity should be near identity"
1310 );
1311
1312 let almost_e = SU2::rotation_x(1e-12);
1313 assert!(
1314 almost_e.is_near_identity(1e-10),
1315 "Small rotation should be near identity"
1316 );
1317
1318 let g = SU2::rotation_y(0.5);
1319 assert!(
1320 !g.is_near_identity(1e-10),
1321 "Large rotation should not be near identity"
1322 );
1323 }
1324
1325 #[test]
1326 fn test_lie_group_generic_algorithm() {
1327 use crate::traits::LieGroup;
1328
1329 // Generic parallel transport function (works for any LieGroup!)
1330 fn parallel_transport<G: LieGroup>(path: &[G]) -> G {
1331 path.iter().fold(G::identity(), |acc, g| acc.compose(g))
1332 }
1333
1334 let path = vec![
1335 SU2::rotation_x(0.1),
1336 SU2::rotation_y(0.2),
1337 SU2::rotation_z(0.3),
1338 ];
1339
1340 let holonomy = parallel_transport(&path);
1341
1342 // Verify it's a valid SU(2) element
1343 assert!(holonomy.verify_unitarity(1e-10));
1344
1345 // Should be non-trivial (not identity)
1346 assert!(holonomy.distance_to_identity() > 0.1);
1347 }
1348
1349 // ========================================================================
1350 // Property-Based Tests for Group Axioms
1351 // ========================================================================
1352 //
1353 // These tests use proptest to verify that SU(2) satisfies the
1354 // mathematical axioms of a Lie group for randomly generated elements.
1355 //
1356 // This is a form of **specification-based testing**: the group axioms
1357 // are the specification, and we verify they hold for all inputs.
1358 //
1359 // Run with: cargo test --features property-tests
1360
1361 #[cfg(feature = "proptest")]
1362 use proptest::prelude::*;
1363
1364 /// Strategy for generating arbitrary SU(2) elements.
1365 ///
1366 /// We generate SU(2) elements by composing three Euler rotations:
1367 /// `g = R_z(α) · R_y(β) · R_x(γ)`
1368 ///
1369 /// This gives a good coverage of SU(2) ≅ S³.
1370 #[cfg(feature = "proptest")]
1371 fn arb_su2() -> impl Strategy<Value = SU2> {
1372 use std::f64::consts::TAU;
1373
1374 // Generate three Euler angles
1375 let alpha = 0.0..TAU;
1376 let beta = 0.0..TAU;
1377 let gamma = 0.0..TAU;
1378
1379 (alpha, beta, gamma).prop_map(|(a, b, c)| {
1380 SU2::rotation_z(a)
1381 .compose(&SU2::rotation_y(b))
1382 .compose(&SU2::rotation_x(c))
1383 })
1384 }
1385
1386 #[cfg(feature = "proptest")]
1387 proptest! {
1388 /// **Group Axiom 1: Identity Element**
1389 ///
1390 /// For all g ∈ SU(2):
1391 /// - e · g = g (left identity)
1392 /// - g · e = g (right identity)
1393 ///
1394 /// where e = identity element
1395 ///
1396 /// Note: We use tolerance 1e-7 to account for floating-point
1397 /// rounding errors in matrix operations.
1398 #[test]
1399 fn prop_identity_axiom(g in arb_su2()) {
1400 let e = SU2::identity();
1401
1402 // Left identity: e · g = g
1403 let left = e.compose(&g);
1404 prop_assert!(
1405 left.distance(&g) < 1e-7,
1406 "Left identity failed: e·g != g, distance = {}",
1407 left.distance(&g)
1408 );
1409
1410 // Right identity: g · e = g
1411 let right = g.compose(&e);
1412 prop_assert!(
1413 right.distance(&g) < 1e-7,
1414 "Right identity failed: g·e != g, distance = {}",
1415 right.distance(&g)
1416 );
1417 }
1418
1419 /// **Group Axiom 2: Inverse Element**
1420 ///
1421 /// For all g ∈ SU(2):
1422 /// - g · g⁻¹ = e (right inverse)
1423 /// - g⁻¹ · g = e (left inverse)
1424 ///
1425 /// where g⁻¹ = inverse of g
1426 ///
1427 /// Note: We use tolerance 1e-7 to account for floating-point
1428 /// rounding errors in matrix operations.
1429 #[test]
1430 fn prop_inverse_axiom(g in arb_su2()) {
1431 let g_inv = g.inverse();
1432
1433 // Right inverse: g · g⁻¹ = e
1434 let right_product = g.compose(&g_inv);
1435 prop_assert!(
1436 right_product.is_near_identity(1e-7),
1437 "Right inverse failed: g·g⁻¹ != e, distance = {}",
1438 right_product.distance_to_identity()
1439 );
1440
1441 // Left inverse: g⁻¹ · g = e
1442 let left_product = g_inv.compose(&g);
1443 prop_assert!(
1444 left_product.is_near_identity(1e-7),
1445 "Left inverse failed: g⁻¹·g != e, distance = {}",
1446 left_product.distance_to_identity()
1447 );
1448 }
1449
1450 /// **Group Axiom 3: Associativity**
1451 ///
1452 /// For all g₁, g₂, g₃ ∈ SU(2):
1453 /// - (g₁ · g₂) · g₃ = g₁ · (g₂ · g₃)
1454 ///
1455 /// Group multiplication is associative.
1456 ///
1457 /// Note: We use tolerance 1e-7 to account for floating-point
1458 /// rounding errors in matrix operations.
1459 #[test]
1460 fn prop_associativity(g1 in arb_su2(), g2 in arb_su2(), g3 in arb_su2()) {
1461 // Left association: (g₁ · g₂) · g₃
1462 let left_assoc = g1.compose(&g2).compose(&g3);
1463
1464 // Right association: g₁ · (g₂ · g₃)
1465 let right_assoc = g1.compose(&g2.compose(&g3));
1466
1467 prop_assert!(
1468 left_assoc.distance(&right_assoc) < 1e-7,
1469 "Associativity failed: (g₁·g₂)·g₃ != g₁·(g₂·g₃), distance = {}",
1470 left_assoc.distance(&right_assoc)
1471 );
1472 }
1473
1474 /// **Lie Group Property: Inverse is Smooth**
1475 ///
1476 /// For SU(2), the inverse operation is smooth (continuously differentiable).
1477 /// We verify this by checking that nearby elements have nearby inverses.
1478 #[test]
1479 fn prop_inverse_continuity(g in arb_su2()) {
1480 // Create a small perturbation
1481 let epsilon = 0.01;
1482 let perturbation = SU2::rotation_x(epsilon);
1483 let g_perturbed = g.compose(&perturbation);
1484
1485 // Check that inverses are close
1486 let inv_distance = g.inverse().distance(&g_perturbed.inverse());
1487
1488 prop_assert!(
1489 inv_distance < 0.1,
1490 "Inverse not continuous: small perturbation caused large inverse change, distance = {}",
1491 inv_distance
1492 );
1493 }
1494
1495 /// **Unitarity Preservation**
1496 ///
1497 /// All SU(2) operations should preserve unitarity.
1498 /// This is not strictly a group axiom, but it's essential for SU(2).
1499 #[test]
1500 fn prop_unitarity_preserved(g1 in arb_su2(), g2 in arb_su2()) {
1501 // Composition preserves unitarity
1502 let product = g1.compose(&g2);
1503 prop_assert!(
1504 product.verify_unitarity(1e-10),
1505 "Composition violated unitarity"
1506 );
1507
1508 // Inverse preserves unitarity
1509 let inv = g1.inverse();
1510 prop_assert!(
1511 inv.verify_unitarity(1e-10),
1512 "Inverse violated unitarity"
1513 );
1514 }
1515
1516 /// **Adjoint Representation: Group Homomorphism**
1517 ///
1518 /// The adjoint representation Ad: G → Aut(𝔤) is a group homomorphism:
1519 /// - Ad_{g₁∘g₂}(X) = Ad_{g₁}(Ad_{g₂}(X))
1520 ///
1521 /// This is a fundamental property that must hold for the adjoint action
1522 /// to be a valid representation of the group.
1523 #[test]
1524 fn prop_adjoint_homomorphism(
1525 g1 in arb_su2(),
1526 g2 in arb_su2(),
1527 x in arb_su2_algebra()
1528 ) {
1529 // Compute Ad_{g₁∘g₂}(X)
1530 let g_composed = g1.compose(&g2);
1531 let left = g_composed.adjoint_action(&x);
1532
1533 // Compute Ad_{g₁}(Ad_{g₂}(X))
1534 let ad_g2_x = g2.adjoint_action(&x);
1535 let right = g1.adjoint_action(&ad_g2_x);
1536
1537 // They should be equal
1538 let diff = left.add(&right.scale(-1.0));
1539 prop_assert!(
1540 diff.norm() < 1e-7,
1541 "Adjoint homomorphism failed: Ad_{{g₁∘g₂}}(X) != Ad_{{g₁}}(Ad_{{g₂}}(X)), diff norm = {}",
1542 diff.norm()
1543 );
1544 }
1545
1546 /// **Adjoint Representation: Identity Action**
1547 ///
1548 /// The identity element acts trivially on the Lie algebra:
1549 /// - Ad_e(X) = X for all X ∈ 𝔤
1550 #[test]
1551 fn prop_adjoint_identity(x in arb_su2_algebra()) {
1552 let e = SU2::identity();
1553 let result = e.adjoint_action(&x);
1554
1555 let diff = result.add(&x.scale(-1.0));
1556 prop_assert!(
1557 diff.norm() < 1e-10,
1558 "Identity action failed: Ad_e(X) != X, diff norm = {}",
1559 diff.norm()
1560 );
1561 }
1562
1563 /// **Adjoint Representation: Lie Bracket Preservation**
1564 ///
1565 /// The adjoint representation preserves the Lie bracket:
1566 /// - Ad_g([X,Y]) = [Ad_g(X), Ad_g(Y)]
1567 ///
1568 /// This is a critical property that ensures the adjoint action
1569 /// is a Lie algebra automorphism.
1570 #[test]
1571 fn prop_adjoint_bracket_preservation(
1572 g in arb_su2(),
1573 x in arb_su2_algebra(),
1574 y in arb_su2_algebra()
1575 ) {
1576 use crate::traits::LieAlgebra;
1577
1578 // Compute Ad_g([X,Y])
1579 let bracket_xy = x.bracket(&y);
1580 let left = g.adjoint_action(&bracket_xy);
1581
1582 // Compute [Ad_g(X), Ad_g(Y)]
1583 let ad_x = g.adjoint_action(&x);
1584 let ad_y = g.adjoint_action(&y);
1585 let right = ad_x.bracket(&ad_y);
1586
1587 // They should be equal
1588 let diff = left.add(&right.scale(-1.0));
1589 prop_assert!(
1590 diff.norm() < 1e-6,
1591 "Bracket preservation failed: Ad_g([X,Y]) != [Ad_g(X), Ad_g(Y)], diff norm = {}",
1592 diff.norm()
1593 );
1594 }
1595
1596 /// **Jacobi Identity: Fundamental Lie Algebra Axiom**
1597 ///
1598 /// The Jacobi identity must hold for all X, Y, Z ∈ 𝔤:
1599 /// - [X, [Y, Z]] + [Y, [Z, X]] + [Z, [X, Y]] = 0
1600 ///
1601 /// This is equivalent to the derivation property of ad_X tested elsewhere,
1602 /// but here we test it directly with three random elements.
1603 #[test]
1604 fn prop_jacobi_identity(
1605 x in arb_su2_algebra(),
1606 y in arb_su2_algebra(),
1607 z in arb_su2_algebra()
1608 ) {
1609 use crate::traits::LieAlgebra;
1610
1611 // Compute [X, [Y, Z]]
1612 let yz = y.bracket(&z);
1613 let term1 = x.bracket(&yz);
1614
1615 // Compute [Y, [Z, X]]
1616 let zx = z.bracket(&x);
1617 let term2 = y.bracket(&zx);
1618
1619 // Compute [Z, [X, Y]]
1620 let xy = x.bracket(&y);
1621 let term3 = z.bracket(&xy);
1622
1623 // Sum should be zero
1624 let sum = term1.add(&term2).add(&term3);
1625 prop_assert!(
1626 sum.norm() < 1e-10,
1627 "Jacobi identity failed: ||[X,[Y,Z]] + [Y,[Z,X]] + [Z,[X,Y]]|| = {:.2e}",
1628 sum.norm()
1629 );
1630 }
1631
1632 /// **Adjoint Representation: Inverse Property**
1633 ///
1634 /// The inverse of an element acts as the inverse transformation:
1635 /// - Ad_{g⁻¹}(Ad_g(X)) = X
1636 #[test]
1637 fn prop_adjoint_inverse(g in arb_su2(), x in arb_su2_algebra()) {
1638 // Apply Ad_g then Ad_{g⁻¹}
1639 let ad_g_x = g.adjoint_action(&x);
1640 let g_inv = g.inverse();
1641 let result = g_inv.adjoint_action(&ad_g_x);
1642
1643 // Should recover X
1644 let diff = result.add(&x.scale(-1.0));
1645 prop_assert!(
1646 diff.norm() < 1e-7,
1647 "Inverse property failed: Ad_{{g⁻¹}}(Ad_g(X)) != X, diff norm = {}",
1648 diff.norm()
1649 );
1650 }
1651 }
1652
1653 /// Strategy for generating arbitrary `Su2Algebra` elements.
1654 ///
1655 /// We generate algebra elements by picking random coefficients in [-π, π]
1656 /// for each of the three basis directions (`σ_x`, `σ_y`, `σ_z`).
1657 #[cfg(feature = "proptest")]
1658 fn arb_su2_algebra() -> impl Strategy<Value = Su2Algebra> {
1659 use proptest::prelude::*;
1660 use std::f64::consts::PI;
1661
1662 ((-PI..PI), (-PI..PI), (-PI..PI)).prop_map(|(a, b, c)| Su2Algebra([a, b, c]))
1663 }
1664
1665 #[test]
1666 #[cfg(feature = "rand")]
1667 fn test_random_haar_unitarity() {
1668 use rand::SeedableRng;
1669
1670 let mut rng = rand::rngs::StdRng::seed_from_u64(42);
1671
1672 // Generate many random elements and verify they're all unitary
1673 for _ in 0..100 {
1674 let g = SU2::random_haar(&mut rng);
1675 assert!(
1676 g.verify_unitarity(1e-10),
1677 "Random Haar element should be unitary"
1678 );
1679 }
1680 }
1681
1682 #[test]
1683 fn test_bracket_bilinearity() {
1684 use crate::traits::LieAlgebra;
1685
1686 let x = Su2Algebra([1.0, 0.0, 0.0]);
1687 let y = Su2Algebra([0.0, 1.0, 0.0]);
1688 let z = Su2Algebra([0.0, 0.0, 1.0]);
1689 let alpha = 2.5;
1690
1691 // Left linearity: [αX + Y, Z] = α[X, Z] + [Y, Z]
1692 let lhs = x.scale(alpha).add(&y).bracket(&z);
1693 let rhs = x.bracket(&z).scale(alpha).add(&y.bracket(&z));
1694 for i in 0..3 {
1695 assert!(
1696 (lhs.0[i] - rhs.0[i]).abs() < 1e-14,
1697 "Left bilinearity failed at component {}: {} vs {}",
1698 i,
1699 lhs.0[i],
1700 rhs.0[i]
1701 );
1702 }
1703
1704 // Right linearity: [Z, αX + Y] = α[Z, X] + [Z, Y]
1705 let lhs = z.bracket(&x.scale(alpha).add(&y));
1706 let rhs = z.bracket(&x).scale(alpha).add(&z.bracket(&y));
1707 for i in 0..3 {
1708 assert!(
1709 (lhs.0[i] - rhs.0[i]).abs() < 1e-14,
1710 "Right bilinearity failed at component {}: {} vs {}",
1711 i,
1712 lhs.0[i],
1713 rhs.0[i]
1714 );
1715 }
1716 }
1717
1718 /// **Jacobi Identity Test with Random Elements**
1719 ///
1720 /// The Jacobi identity is the defining axiom of a Lie algebra:
1721 /// ```text
1722 /// [X, [Y, Z]] + [Y, [Z, X]] + [Z, [X, Y]] = 0
1723 /// ```
1724 ///
1725 /// This test verifies it holds for random algebra elements.
1726 #[test]
1727 fn test_jacobi_identity_random() {
1728 use crate::traits::LieAlgebra;
1729 use rand::SeedableRng;
1730 use rand_distr::{Distribution, StandardNormal};
1731
1732 let mut rng = rand::rngs::StdRng::seed_from_u64(12345);
1733
1734 // Test with 100 random triples
1735 for _ in 0..100 {
1736 // Generate random algebra elements
1737 let x = Su2Algebra([
1738 StandardNormal.sample(&mut rng),
1739 StandardNormal.sample(&mut rng),
1740 StandardNormal.sample(&mut rng),
1741 ]);
1742 let y = Su2Algebra([
1743 StandardNormal.sample(&mut rng),
1744 StandardNormal.sample(&mut rng),
1745 StandardNormal.sample(&mut rng),
1746 ]);
1747 let z = Su2Algebra([
1748 StandardNormal.sample(&mut rng),
1749 StandardNormal.sample(&mut rng),
1750 StandardNormal.sample(&mut rng),
1751 ]);
1752
1753 // Compute [X, [Y, Z]]
1754 let yz = y.bracket(&z);
1755 let term1 = x.bracket(&yz);
1756
1757 // Compute [Y, [Z, X]]
1758 let zx = z.bracket(&x);
1759 let term2 = y.bracket(&zx);
1760
1761 // Compute [Z, [X, Y]]
1762 let xy = x.bracket(&y);
1763 let term3 = z.bracket(&xy);
1764
1765 // Sum should be zero
1766 let sum = term1.add(&term2).add(&term3);
1767 assert!(
1768 sum.norm() < 1e-10,
1769 "Jacobi identity failed: ||[X,[Y,Z]] + [Y,[Z,X]] + [Z,[X,Y]]|| = {:.2e}",
1770 sum.norm()
1771 );
1772 }
1773 }
1774
1775 /// **Exp-Log Round-Trip Test**
1776 ///
1777 /// For any algebra element X with ||X|| < π, we should have:
1778 /// ```text
1779 /// log(exp(X)) = X
1780 /// ```
1781 #[test]
1782 fn test_exp_log_roundtrip() {
1783 use crate::traits::{LieAlgebra, LieGroup};
1784 use rand::SeedableRng;
1785 use rand_distr::{Distribution, Uniform};
1786
1787 let mut rng = rand::rngs::StdRng::seed_from_u64(54321);
1788 let dist = Uniform::new(-2.0, 2.0); // Stay within log domain
1789
1790 for _ in 0..100 {
1791 let x = Su2Algebra([
1792 dist.sample(&mut rng),
1793 dist.sample(&mut rng),
1794 dist.sample(&mut rng),
1795 ]);
1796
1797 // exp then log
1798 let g = SU2::exp(&x);
1799 let x_recovered = g.log().expect("log should succeed for exp output");
1800
1801 // Should recover original (up to numerical precision)
1802 let diff = x.add(&x_recovered.scale(-1.0));
1803 assert!(
1804 diff.norm() < 1e-10,
1805 "log(exp(X)) should equal X: ||diff|| = {:.2e}",
1806 diff.norm()
1807 );
1808 }
1809 }
1810
1811 /// **Log-Exp Round-Trip Test**
1812 ///
1813 /// For any group element g, we should have:
1814 /// ```text
1815 /// exp(log(g)) = g
1816 /// ```
1817 #[test]
1818 #[cfg(feature = "rand")]
1819 fn test_log_exp_roundtrip() {
1820 use crate::traits::LieGroup;
1821 use rand::SeedableRng;
1822
1823 let mut rng = rand::rngs::StdRng::seed_from_u64(98765);
1824
1825 for _ in 0..100 {
1826 let g = SU2::random_haar(&mut rng);
1827
1828 // log then exp
1829 let x = g.log().expect("log should succeed for valid SU(2) element");
1830 let g_recovered = SU2::exp(&x);
1831
1832 // Should recover original (numerical precision varies with rotation angle)
1833 let diff = g.compose(&g_recovered.inverse()).distance_to_identity();
1834 assert!(
1835 diff < 1e-7,
1836 "exp(log(g)) should equal g: diff = {:.2e}",
1837 diff
1838 );
1839 }
1840 }
1841
1842 #[test]
1843 #[cfg(feature = "rand")]
1844 fn test_random_haar_distribution() {
1845 use rand::SeedableRng;
1846
1847 let mut rng = rand::rngs::StdRng::seed_from_u64(123);
1848
1849 // Generate many samples and check statistical properties
1850 let samples: Vec<SU2> = (0..1000).map(|_| SU2::random_haar(&mut rng)).collect();
1851
1852 // Check that average distance to identity is approximately correct
1853 // For uniform distribution on SU(2) ≅ S³, the expected distance is non-trivial
1854 let mean_distance: f64 =
1855 samples.iter().map(SU2::distance_to_identity).sum::<f64>() / samples.len() as f64;
1856
1857 // For S³ with uniform (Haar) measure, the mean geodesic distance from identity
1858 // is approximately π (the maximum is π for antipodal points).
1859 // Empirically, the mean is close to π.
1860 assert!(
1861 mean_distance > 2.5 && mean_distance < 3.5,
1862 "Mean distance from identity should be ~π, got {}",
1863 mean_distance
1864 );
1865
1866 // Check that some elements are far from identity
1867 let far_from_identity = samples
1868 .iter()
1869 .filter(|g| g.distance_to_identity() > std::f64::consts::PI / 2.0)
1870 .count();
1871
1872 assert!(
1873 far_from_identity > 100,
1874 "Should have many elements far from identity, got {}",
1875 far_from_identity
1876 );
1877 }
1878
1879 #[test]
1880 fn test_adjoint_action_simple() {
1881 use crate::traits::LieGroup;
1882
1883 // Test with identity: Ad_e(X) = X
1884 let e = SU2::identity();
1885 let x = Su2Algebra([1.0, 0.0, 0.0]);
1886 let result = e.adjoint_action(&x);
1887
1888 println!("Identity test: X = {:?}, Ad_e(X) = {:?}", x, result);
1889 assert!((result.0[0] - x.0[0]).abs() < 1e-10);
1890 assert!((result.0[1] - x.0[1]).abs() < 1e-10);
1891 assert!((result.0[2] - x.0[2]).abs() < 1e-10);
1892
1893 // Test with rotation around Z by 90 degrees
1894 // Should rotate X basis element into Y basis element
1895 let g = SU2::rotation_z(std::f64::consts::FRAC_PI_2);
1896 let x_basis = Su2Algebra([1.0, 0.0, 0.0]);
1897 let rotated = g.adjoint_action(&x_basis);
1898
1899 println!(
1900 "Rotation test: X = {:?}, Ad_{{Rz(π/2)}}(X) = {:?}",
1901 x_basis, rotated
1902 );
1903 // Should be approximately (0, 1, 0) - rotated 90° around Z
1904 }
1905
1906 // ========================================================================
1907 // Casimir Operator Tests
1908 // ========================================================================
1909
1910 #[test]
1911 fn test_su2_casimir_scalar() {
1912 // j = 0 (scalar): c₂ = 0
1913 use crate::representation::Spin;
1914 use crate::Casimir;
1915
1916 let c2 = Su2Algebra::quadratic_casimir_eigenvalue(&Spin::ZERO);
1917 assert_eq!(c2, 0.0, "Casimir of scalar representation should be 0");
1918 }
1919
1920 #[test]
1921 fn test_su2_casimir_spinor() {
1922 // j = 1/2 (spinor): c₂ = 3/4
1923 use crate::representation::Spin;
1924 use crate::Casimir;
1925
1926 let c2 = Su2Algebra::quadratic_casimir_eigenvalue(&Spin::HALF);
1927 assert_eq!(c2, 0.75, "Casimir of spinor representation should be 3/4");
1928 }
1929
1930 #[test]
1931 fn test_su2_casimir_vector() {
1932 // j = 1 (vector/adjoint): c₂ = 2
1933 use crate::representation::Spin;
1934 use crate::Casimir;
1935
1936 let c2 = Su2Algebra::quadratic_casimir_eigenvalue(&Spin::ONE);
1937 assert_eq!(c2, 2.0, "Casimir of vector representation should be 2");
1938 }
1939
1940 #[test]
1941 fn test_su2_casimir_j_three_halves() {
1942 // j = 3/2: c₂ = 15/4
1943 use crate::representation::Spin;
1944 use crate::Casimir;
1945
1946 let j_three_halves = Spin::from_half_integer(3);
1947 let c2 = Su2Algebra::quadratic_casimir_eigenvalue(&j_three_halves);
1948 assert_eq!(c2, 3.75, "Casimir for j=3/2 should be 15/4 = 3.75");
1949 }
1950
1951 #[test]
1952 fn test_su2_casimir_formula() {
1953 // Test the general formula c₂(j) = j(j+1) for various spins
1954 use crate::representation::Spin;
1955 use crate::Casimir;
1956
1957 for two_j in 0..10 {
1958 let spin = Spin::from_half_integer(two_j);
1959 let j = spin.value();
1960 let c2 = Su2Algebra::quadratic_casimir_eigenvalue(&spin);
1961 let expected = j * (j + 1.0);
1962
1963 assert_relative_eq!(c2, expected, epsilon = 1e-10);
1964 }
1965 }
1966
1967 #[test]
1968 fn test_su2_rank() {
1969 // SU(2) has rank 1
1970 use crate::Casimir;
1971
1972 assert_eq!(Su2Algebra::rank(), 1, "SU(2) should have rank 1");
1973 assert_eq!(
1974 Su2Algebra::num_casimirs(),
1975 1,
1976 "SU(2) should have 1 Casimir operator"
1977 );
1978 }
1979
1980 // ========================================================================
1981 // Stable Logarithm Tests
1982 // ========================================================================
1983
1984 #[test]
1985 fn test_log_stable_identity() {
1986 let e = SU2::identity();
1987 let log_e = e.log_stable().expect("log of identity should succeed");
1988 assert!(log_e.norm() < 1e-10, "log(I) should be zero");
1989 }
1990
1991 #[test]
1992 fn test_log_stable_small_rotation() {
1993 let theta = 0.1;
1994 let g = SU2::rotation_x(theta);
1995 let log_g = g
1996 .log_stable()
1997 .expect("log of small rotation should succeed");
1998
1999 // Should give θ * (1, 0, 0)
2000 assert!((log_g.0[0] - theta).abs() < 1e-10);
2001 assert!(log_g.0[1].abs() < 1e-10);
2002 assert!(log_g.0[2].abs() < 1e-10);
2003 }
2004
2005 #[test]
2006 fn test_log_stable_large_rotation() {
2007 let theta = 2.5; // Less than π
2008 let g = SU2::rotation_y(theta);
2009 let log_g = g.log_stable().expect("log of rotation < π should succeed");
2010
2011 // Should give θ * (0, 1, 0)
2012 assert!(log_g.0[0].abs() < 1e-10);
2013 assert!((log_g.0[1] - theta).abs() < 1e-10);
2014 assert!(log_g.0[2].abs() < 1e-10);
2015 }
2016
2017 #[test]
2018 fn test_log_stable_vs_log_consistency() {
2019 use crate::traits::{LieAlgebra, LieGroup};
2020
2021 // For rotations away from 2π, log_stable and log should agree
2022 // Note: θ = π is NOT a singularity for SU(2), so we include it.
2023 for theta in [0.1, 0.5, 1.0, 1.5, 2.0, 2.5, std::f64::consts::PI] {
2024 let g = SU2::rotation_z(theta);
2025 let log_standard = g.log().expect("log should succeed");
2026 let log_stable = g.log_stable().expect("log_stable should succeed");
2027
2028 let diff = log_standard.add(&log_stable.scale(-1.0)).norm();
2029 assert!(
2030 diff < 1e-10,
2031 "log and log_stable should agree for θ = {}: diff = {:.2e}",
2032 theta,
2033 diff
2034 );
2035 }
2036 }
2037
2038 #[test]
2039 fn test_log_with_condition_returns_condition() {
2040 let theta = 1.0;
2041 let g = SU2::rotation_x(theta);
2042 let (log_g, cond) = g
2043 .log_with_condition()
2044 .expect("log_with_condition should succeed");
2045
2046 // Should be well-conditioned for θ = 1
2047 assert!(
2048 cond.is_well_conditioned(),
2049 "θ = 1 should be well-conditioned"
2050 );
2051 assert!(
2052 (cond.angle() - theta).abs() < 1e-10,
2053 "reported angle should match"
2054 );
2055 assert!((log_g.0[0] - theta).abs() < 1e-10);
2056 }
2057
2058 #[test]
2059 fn test_log_with_condition_near_pi() {
2060 // At θ ≈ π, sin(θ/2) ≈ 1, so axis extraction is perfectly stable.
2061 // The cut locus for SU(2) is at θ = 2π (U = -I), not θ = π.
2062 let theta = std::f64::consts::PI - 0.01;
2063 let g = SU2::rotation_z(theta);
2064 let (log_g, cond) = g
2065 .log_with_condition()
2066 .expect("log_with_condition should return best-effort");
2067
2068 // At θ ≈ π, sin(θ/2) ≈ 1, so numerical extraction is stable
2069 assert!(
2070 cond.is_well_conditioned(),
2071 "θ ≈ π should be numerically well-conditioned: κ = {}",
2072 cond.condition_number()
2073 );
2074
2075 // Distance to cut locus (2π) should be about π
2076 assert!(
2077 (cond.distance_to_cut_locus() - (std::f64::consts::PI + 0.01)).abs() < 1e-10,
2078 "distance to cut locus should be ≈ π + 0.01: got {}",
2079 cond.distance_to_cut_locus()
2080 );
2081
2082 // Result should be correct
2083 assert!((log_g.0[2] - theta).abs() < 1e-6);
2084 }
2085
2086 #[test]
2087 fn test_log_condition_from_angle() {
2088 use crate::{LogCondition, LogQuality};
2089
2090 // Small angle θ = 0.5: sin(0.25) ≈ 0.247, κ ≈ 4.0 → Good
2091 // The condition number tracks stability of axis extraction (dividing by sin(θ/2))
2092 let cond_small = LogCondition::from_angle(0.5);
2093 assert_eq!(cond_small.quality(), LogQuality::Good);
2094 assert!(cond_small.condition_number() > 2.0 && cond_small.condition_number() < 10.0);
2095
2096 // π/2: sin(π/4) ≈ 0.707, κ ≈ 1.4 → Excellent
2097 let cond_half_pi = LogCondition::from_angle(std::f64::consts::FRAC_PI_2);
2098 assert_eq!(cond_half_pi.quality(), LogQuality::Excellent);
2099 assert!(cond_half_pi.is_well_conditioned());
2100
2101 // Near π: sin((π-0.001)/2) ≈ 1, κ ≈ 1.0 → Excellent (numerically stable)
2102 // Note: This is NUMERICALLY stable even though the axis is MATHEMATICALLY ambiguous
2103 let cond_near_pi = LogCondition::from_angle(std::f64::consts::PI - 0.001);
2104 assert!(
2105 cond_near_pi.is_well_conditioned(),
2106 "Near π should be numerically stable: κ = {}",
2107 cond_near_pi.condition_number()
2108 );
2109
2110 // Near zero: sin(0.001/2) ≈ 0.0005, κ ≈ 2000 → AtSingularity
2111 // This is where axis extraction becomes unstable (dividing by small number)
2112 let cond_near_zero = LogCondition::from_angle(0.001);
2113 assert!(
2114 !cond_near_zero.is_well_conditioned(),
2115 "Near zero should have poor conditioning: κ = {}",
2116 cond_near_zero.condition_number()
2117 );
2118 }
2119
2120 #[test]
2121 fn test_log_stable_roundtrip() {
2122 use crate::traits::{LieAlgebra, LieGroup};
2123 use rand::SeedableRng;
2124 use rand_distr::{Distribution, Uniform};
2125
2126 let mut rng = rand::rngs::StdRng::seed_from_u64(99999);
2127 let dist = Uniform::new(-2.5, 2.5); // Stay within log domain but include large angles
2128
2129 for _ in 0..100 {
2130 let x = Su2Algebra([
2131 dist.sample(&mut rng),
2132 dist.sample(&mut rng),
2133 dist.sample(&mut rng),
2134 ]);
2135
2136 // Skip elements near the singularity
2137 if x.norm() > std::f64::consts::PI - 0.2 {
2138 continue;
2139 }
2140
2141 // exp then log_stable
2142 let g = SU2::exp(&x);
2143 let x_recovered = g.log_stable().expect("log_stable should succeed");
2144
2145 // Should recover original
2146 let diff = x.add(&x_recovered.scale(-1.0));
2147 assert!(
2148 diff.norm() < 1e-8,
2149 "log_stable(exp(X)) should equal X: ||diff|| = {:.2e}",
2150 diff.norm()
2151 );
2152 }
2153 }
2154
2155 #[test]
2156 fn test_log_quality_display() {
2157 use crate::LogQuality;
2158
2159 assert_eq!(format!("{}", LogQuality::Excellent), "Excellent");
2160 assert_eq!(format!("{}", LogQuality::Good), "Good");
2161 assert_eq!(format!("{}", LogQuality::Acceptable), "Acceptable");
2162 assert_eq!(format!("{}", LogQuality::Poor), "Poor");
2163 assert_eq!(format!("{}", LogQuality::AtSingularity), "AtSingularity");
2164 }
2165}