Skip to main content

lie_groups/
quaternion.rs

1//! Quaternionic Formulation of SU(2)
2//!
3//! **Addressing Penrose's Recommendation:**
4//! > "SU(2) ≅ unit quaternions. Why not use this beautiful structure?"
5//!
6//! This module implements SU(2) using quaternions instead of 2×2 matrices.
7//!
8//! # Mathematical Background
9//!
10//! **Quaternions:** ℍ = {a + bi + cj + dk | a,b,c,d ∈ ℝ}
11//! with multiplication rules: i² = j² = k² = ijk = -1
12//!
13//! **Unit Quaternions:** S³ = {q ∈ ℍ | |q| = 1} ⊂ ℍ
14//!
15//! **Fundamental Isomorphism:** SU(2) ≅ S³
16//!
17//! The map is:
18//! ```text
19//! SU(2) matrix [[α, -β*], [β, α*]]  ↔  Quaternion q = α + βj
20//!                                        where α = a + bi, β = c + di
21//!                                        ⟺ q = a + bi + cj + dk
22//! ```
23//!
24//! # Advantages Over Matrix Representation
25//!
26//! 1. **Compact**: 4 reals vs 8 reals (complex 2×2 matrix)
27//! 2. **Efficient**: Quaternion multiplication is faster than matrix multiplication
28//! 3. **Geometric**: Direct axis-angle interpretation for rotations
29//! 4. **Numerical**: Better numerical stability (no matrix decomposition needed)
30//! 5. **Elegant**: Exponential map is cleaner
31//!
32//! # The SU(2) Double Cover of SO(3) ⭐
33//!
34//! **Critical Physics Insight:**
35//!
36//! SU(2) is a **double cover** of SO(3) (3D rotations):
37//! ```text
38//! SU(2) → SO(3)  is a 2-to-1 surjective homomorphism
39//! ```
40//!
41//! **Meaning:** Quaternions **q** and **-q** represent the **same rotation** in 3D space,
42//! but **different quantum states** in SU(2).
43//!
44//! ## Mathematical Statement
45//!
46//! For any unit quaternion q:
47//! - **q** rotates a vector v by angle θ around axis n̂
48//! - **-q** also rotates v by angle θ around axis n̂  (same rotation!)
49//! - But q ≠ -q as group elements in SU(2)
50//!
51//! Path from q → -q → q (going 4π around) is **topologically non-trivial**.
52//!
53//! ## Physical Consequences
54//!
55//! 1. **Fermions (spin-½ particles)**:
56//!    - Under 2π rotation: ψ → -ψ  (sign flip!)
57//!    - Under 4π rotation: ψ → ψ  (back to original)
58//!    - This is why fermions have **half-integer spin**
59//!
60//! 2. **Bosons (integer spin)**:
61//!    - Under 2π rotation: ψ → ψ  (no sign change)
62//!    - Described by SO(3), not SU(2)
63//!
64//! 3. **Berry Phase**:
65//!    - Geometric phase accumulated when traversing closed path in SU(2)
66//!    - q → -q → q acquires phase π (Pancharatnam-Berry phase)
67//!
68//! ## Example: Rotation by 2π
69//!
70//! ```text
71//! Rotation by angle θ = 2π around z-axis:
72//!   q = cos(π) + sin(π)·k = -1 + 0·k = -1  (minus identity!)
73//!
74//! But -1 acts on vector v as:
75//!   v → -1·v·(-1)⁻¹ = -1·v·(-1) = v  (identity rotation in SO(3))
76//! ```
77//!
78//! ## Gauge Theory Interpretation
79//!
80//! In gauge theory on networks:
81//! - **Connections** A_{ij} ∈ SU(2) carry full SU(2) structure (not just rotations)
82//! - **Wilson loops** W(C) can wind around non-trivially in SU(2)
83//! - Path C → 2π rotation can give W(C) = -I ≠ I (topological effect!)
84//! - This is invisible in classical physics but crucial for quantum systems
85//!
86//! ## References
87//!
88//! - Dirac, P.A.M.: "The Principles of Quantum Mechanics" (1930) - Original treatment
89//! - Feynman & Hibbs: "Quantum Mechanics and Path Integrals" - Spin-statistics connection
90//! - Penrose: "The Road to Reality" § 11.3 - Geometric interpretation
91//! - Aharonov & Susskind: "Observability of the sign change of spinors" (1967)
92
93use std::ops::{Mul, MulAssign};
94
95/// Unit quaternion representing an element of SU(2)
96///
97/// Quaternion q = w + xi + yj + zk where w² + x² + y² + z² = 1
98///
99/// # Geometric Interpretation
100///
101/// Every unit quaternion represents a rotation in 3D space:
102/// - Rotation by angle θ around axis n̂ = (nx, ny, nz):
103///   q = cos(θ/2) + sin(θ/2)(nx·i + ny·j + nz·k)
104///
105/// # ⚠️ Important: Double Cover Property
106///
107/// **q and -q represent the SAME rotation but DIFFERENT quantum states:**
108/// - Both rotate 3D vectors identically (v → q·v·q⁻¹ = (-q)·v·(-q)⁻¹)
109/// - But q ≠ -q as SU(2) group elements
110/// - Consequence: Fermions change sign under 2π rotation (ψ → -ψ)
111///
112/// This is NOT a numerical issue or degeneracy - it's fundamental topology!
113///
114/// # Relation to SU(2) Matrix
115///
116/// q = w + xi + yj + zk corresponds to SU(2) matrix:
117/// ```text
118/// U = [[ w + ix,  -y + iz],
119///      [ y + iz,   w - ix]]
120/// ```
121///
122/// Note: -q gives matrix -U, which acts identically on vectors but represents
123/// a different element of SU(2) (relevant for spinors/fermions).
124#[derive(Debug, Clone, Copy, PartialEq)]
125pub struct UnitQuaternion {
126    // Private fields enforce unit norm invariant through constructors
127    w: f64,
128    x: f64,
129    y: f64,
130    z: f64,
131}
132
133impl UnitQuaternion {
134    /// Get the real (scalar) component.
135    #[inline]
136    #[must_use]
137    pub fn w(&self) -> f64 {
138        self.w
139    }
140
141    /// Get the imaginary i component.
142    #[inline]
143    #[must_use]
144    pub fn x(&self) -> f64 {
145        self.x
146    }
147
148    /// Get the imaginary j component.
149    #[inline]
150    #[must_use]
151    pub fn y(&self) -> f64 {
152        self.y
153    }
154
155    /// Get the imaginary k component.
156    #[inline]
157    #[must_use]
158    pub fn z(&self) -> f64 {
159        self.z
160    }
161
162    /// Get all components as a tuple (w, x, y, z).
163    #[inline]
164    #[must_use]
165    pub fn components(&self) -> (f64, f64, f64, f64) {
166        (self.w, self.x, self.y, self.z)
167    }
168
169    /// Get all components as an array [w, x, y, z].
170    #[inline]
171    #[must_use]
172    pub fn to_array(&self) -> [f64; 4] {
173        [self.w, self.x, self.y, self.z]
174    }
175}
176
177// ============================================================================
178// Numerical Tolerance Constants
179// ============================================================================
180//
181// These thresholds are chosen based on IEEE 754 f64 precision (ε ≈ 2.2e-16)
182// and practical considerations for accumulated floating-point error.
183//
184// Rationale: 1e-10 = ~10^6 × machine_epsilon
185// This allows for approximately 6 orders of magnitude of accumulated error
186// from typical operations (multiplication, addition, normalization) before
187// triggering degeneracy handling.
188//
189// For quaternion operations on unit quaternions, typical error accumulation:
190// - Single operation: O(ε) ≈ 2e-16
191// - After ~1000 operations: O(1000ε) ≈ 2e-13
192// - Threshold 1e-10 provides ~1000× safety margin
193//
194// Reference: Higham, "Accuracy and Stability of Numerical Algorithms" (2002)
195
196/// Threshold below which a norm is considered zero (degenerate case).
197/// Used for normalization guards and axis-angle singularity detection.
198const NORM_EPSILON: f64 = 1e-10;
199
200/// Threshold for detecting near-identity quaternions.
201/// Used in logarithm to avoid division by near-zero sin(θ/2).
202const IDENTITY_EPSILON: f64 = 1e-10;
203
204/// Threshold for detecting rotation angle near π (axis ambiguity).
205/// At θ = π, the rotation axis becomes undefined (any axis works).
206/// Reserved for future implementation of θ≈π handling.
207#[allow(dead_code)]
208const SINGULARITY_EPSILON: f64 = 1e-6;
209
210impl UnitQuaternion {
211    /// Identity element: q = 1
212    #[must_use]
213    pub fn identity() -> Self {
214        Self {
215            w: 1.0,
216            x: 0.0,
217            y: 0.0,
218            z: 0.0,
219        }
220    }
221
222    /// Create from components (automatically normalizes to unit quaternion)
223    ///
224    /// Returns identity if input norm < `NORM_EPSILON` (degenerate case).
225    #[must_use]
226    pub fn new(w: f64, x: f64, y: f64, z: f64) -> Self {
227        let norm = (w * w + x * x + y * y + z * z).sqrt();
228        if norm < NORM_EPSILON {
229            return Self::identity();
230        }
231        Self {
232            w: w / norm,
233            x: x / norm,
234            y: y / norm,
235            z: z / norm,
236        }
237    }
238
239    /// Create from axis-angle representation
240    ///
241    /// # Arguments
242    /// * `axis` - Rotation axis (nx, ny, nz), will be normalized
243    /// * `angle` - Rotation angle in radians
244    ///
245    /// # Returns
246    /// Unit quaternion q = cos(θ/2) + sin(θ/2)(nx·i + ny·j + nz·k)
247    #[must_use]
248    pub fn from_axis_angle(axis: [f64; 3], angle: f64) -> Self {
249        let [nx, ny, nz] = axis;
250        let axis_norm = (nx * nx + ny * ny + nz * nz).sqrt();
251
252        // Degenerate axis: return identity (no rotation)
253        if axis_norm < NORM_EPSILON {
254            return Self::identity();
255        }
256
257        let half_angle = angle / 2.0;
258        let (sin_half, cos_half) = half_angle.sin_cos();
259
260        Self {
261            w: cos_half,
262            x: sin_half * nx / axis_norm,
263            y: sin_half * ny / axis_norm,
264            z: sin_half * nz / axis_norm,
265        }
266    }
267
268    /// Extract axis and angle from quaternion
269    ///
270    /// # Returns
271    /// (axis, angle) where axis is normalized and angle ∈ [0, 2π]
272    #[must_use]
273    pub fn to_axis_angle(&self) -> ([f64; 3], f64) {
274        // q = cos(θ/2) + sin(θ/2)·n̂
275        let angle = 2.0 * self.w.clamp(-1.0, 1.0).acos();
276        let sin_half = (1.0 - self.w * self.w).max(0.0).sqrt();
277
278        if sin_half < IDENTITY_EPSILON {
279            // Near identity, axis is arbitrary
280            return ([1.0, 0.0, 0.0], 0.0);
281        }
282
283        let axis = [self.x / sin_half, self.y / sin_half, self.z / sin_half];
284
285        (axis, angle)
286    }
287
288    /// Rotation around X-axis by angle θ
289    #[must_use]
290    pub fn rotation_x(theta: f64) -> Self {
291        Self::from_axis_angle([1.0, 0.0, 0.0], theta)
292    }
293
294    /// Rotation around Y-axis by angle θ
295    #[must_use]
296    pub fn rotation_y(theta: f64) -> Self {
297        Self::from_axis_angle([0.0, 1.0, 0.0], theta)
298    }
299
300    /// Rotation around Z-axis by angle θ
301    #[must_use]
302    pub fn rotation_z(theta: f64) -> Self {
303        Self::from_axis_angle([0.0, 0.0, 1.0], theta)
304    }
305
306    /// Quaternion conjugate: q* = w - xi - yj - zk
307    ///
308    /// For unit quaternions: q* = q⁻¹ (inverse)
309    #[must_use]
310    pub fn conjugate(&self) -> Self {
311        Self {
312            w: self.w,
313            x: -self.x,
314            y: -self.y,
315            z: -self.z,
316        }
317    }
318
319    /// Quaternion inverse: q⁻¹
320    ///
321    /// For unit quaternions: q⁻¹ = q*
322    #[must_use]
323    pub fn inverse(&self) -> Self {
324        self.conjugate()
325    }
326
327    /// Norm squared: |q|² = w² + x² + y² + z²
328    ///
329    /// Should always be 1 for unit quaternions
330    #[must_use]
331    pub fn norm_squared(&self) -> f64 {
332        self.w * self.w + self.x * self.x + self.y * self.y + self.z * self.z
333    }
334
335    /// Norm: |q| = sqrt(w² + x² + y² + z²)
336    ///
337    /// Should always be 1 for unit quaternions
338    #[must_use]
339    pub fn norm(&self) -> f64 {
340        self.norm_squared().sqrt()
341    }
342
343    /// Renormalize to unit quaternion (fix numerical drift)
344    #[must_use]
345    pub fn normalize(&self) -> Self {
346        Self::new(self.w, self.x, self.y, self.z)
347    }
348
349    /// Distance to identity (geodesic distance on S³)
350    ///
351    /// d(q, 1) = 2·arccos(|w|)
352    ///
353    /// This is the rotation angle in [0, π]
354    #[must_use]
355    pub fn distance_to_identity(&self) -> f64 {
356        2.0 * self.w.abs().clamp(-1.0, 1.0).acos()
357    }
358
359    /// Geodesic distance to another quaternion
360    ///
361    /// d(q₁, q₂) = arccos(|q₁·q₂|) where · is quaternion dot product
362    #[must_use]
363    pub fn distance_to(&self, other: &Self) -> f64 {
364        let dot = self.w * other.w + self.x * other.x + self.y * other.y + self.z * other.z;
365        dot.abs().clamp(-1.0, 1.0).acos()
366    }
367
368    /// Spherical linear interpolation (SLERP)
369    ///
370    /// Interpolate between self and other with parameter t ∈ `[0,1]`
371    /// Returns shortest path on S³
372    #[must_use]
373    pub fn slerp(&self, other: &Self, t: f64) -> Self {
374        let dot = self.w * other.w + self.x * other.x + self.y * other.y + self.z * other.z;
375
376        // If quaternions are close, use linear interpolation
377        if dot.abs() > 0.9995 {
378            return Self::new(
379                self.w + t * (other.w - self.w),
380                self.x + t * (other.x - self.x),
381                self.y + t * (other.y - self.y),
382                self.z + t * (other.z - self.z),
383            );
384        }
385
386        // Ensure shortest path
387        let (other_w, other_x, other_y, other_z) = if dot < 0.0 {
388            (-other.w, -other.x, -other.y, -other.z)
389        } else {
390            (other.w, other.x, other.y, other.z)
391        };
392
393        let theta = dot.abs().acos();
394        let sin_theta = theta.sin();
395
396        let a = ((1.0 - t) * theta).sin() / sin_theta;
397        let b = (t * theta).sin() / sin_theta;
398
399        Self::new(
400            a * self.w + b * other_w,
401            a * self.x + b * other_x,
402            a * self.y + b * other_y,
403            a * self.z + b * other_z,
404        )
405    }
406
407    /// Exponential map from Lie algebra 𝔰𝔲(2) to group SU(2)
408    ///
409    /// Given a vector v = (v₁, v₂, v₃) ∈ ℝ³ (Lie algebra element),
410    /// compute exp(v) as a unit quaternion.
411    ///
412    /// Formula: exp(v) = cos(|v|/2) + sin(|v|/2)·(v/|v|)
413    #[must_use]
414    pub fn exp(v: [f64; 3]) -> Self {
415        let norm = (v[0] * v[0] + v[1] * v[1] + v[2] * v[2]).sqrt();
416
417        if norm < NORM_EPSILON {
418            return Self::identity();
419        }
420
421        let half_norm = norm / 2.0;
422        let (sin_half, cos_half) = half_norm.sin_cos();
423
424        Self {
425            w: cos_half,
426            x: sin_half * v[0] / norm,
427            y: sin_half * v[1] / norm,
428            z: sin_half * v[2] / norm,
429        }
430    }
431
432    /// Logarithm map from group SU(2) to Lie algebra 𝔰𝔲(2)
433    ///
434    /// Returns vector v ∈ ℝ³ such that exp(v) = q
435    ///
436    /// Formula: log(q) = (θ/sin(θ/2))·(x, y, z) where θ = 2·arccos(w)
437    #[must_use]
438    pub fn log(&self) -> [f64; 3] {
439        let theta = 2.0 * self.w.clamp(-1.0, 1.0).acos();
440        let sin_half = (1.0 - self.w * self.w).max(0.0).sqrt();
441
442        if sin_half < IDENTITY_EPSILON {
443            // Near identity
444            return [0.0, 0.0, 0.0];
445        }
446
447        let scale = theta / sin_half;
448        [scale * self.x, scale * self.y, scale * self.z]
449    }
450
451    /// Convert to SU(2) matrix representation (for compatibility)
452    ///
453    /// Returns 2×2 complex matrix:
454    /// U = [[ w + ix,  -y + iz],
455    ///      [ y + iz,   w - ix]]
456    #[must_use]
457    pub fn to_matrix(&self) -> [[num_complex::Complex64; 2]; 2] {
458        use num_complex::Complex64;
459
460        [
461            [
462                Complex64::new(self.w, self.x),
463                Complex64::new(-self.y, self.z),
464            ],
465            [
466                Complex64::new(self.y, self.z),
467                Complex64::new(self.w, -self.x),
468            ],
469        ]
470    }
471
472    /// Create from SU(2) matrix
473    ///
474    /// Expects matrix [[α, -β*], [β, α*]] with |α|² + |β|² = 1
475    #[must_use]
476    pub fn from_matrix(matrix: [[num_complex::Complex64; 2]; 2]) -> Self {
477        // Extract α = a + ib from matrix[0][0]
478        let a = matrix[0][0].re;
479        let b = matrix[0][0].im;
480
481        // Extract β = c + id from matrix[1][0]
482        let c = matrix[1][0].re;
483        let d = matrix[1][0].im;
484
485        // q = a + bi + cj + dk
486        Self::new(a, b, c, d)
487    }
488
489    /// Act on a 3D vector by conjugation: v' = qvq*
490    ///
491    /// This is the rotation action of SU(2) on ℝ³
492    /// Identifies v = (x,y,z) with quaternion xi + yj + zk
493    ///
494    /// Uses the direct formula (Rodrigues) for efficiency:
495    /// v' = v + 2w(n×v) + 2(n×(n×v)) where n = (x,y,z)
496    #[must_use]
497    pub fn rotate_vector(&self, v: [f64; 3]) -> [f64; 3] {
498        let (w, qx, qy, qz) = (self.w, self.x, self.y, self.z);
499        let (vx, vy, vz) = (v[0], v[1], v[2]);
500
501        // Cross product: n × v
502        let cx = qy * vz - qz * vy;
503        let cy = qz * vx - qx * vz;
504        let cz = qx * vy - qy * vx;
505
506        // Cross product: n × (n × v)
507        let ccx = qy * cz - qz * cy;
508        let ccy = qz * cx - qx * cz;
509        let ccz = qx * cy - qy * cx;
510
511        // v' = v + 2w(n×v) + 2(n×(n×v))
512        [
513            vx + 2.0 * (w * cx + ccx),
514            vy + 2.0 * (w * cy + ccy),
515            vz + 2.0 * (w * cz + ccz),
516        ]
517    }
518
519    /// Verify this is approximately a unit quaternion
520    #[must_use]
521    pub fn verify_unit(&self, tolerance: f64) -> bool {
522        (self.norm_squared() - 1.0).abs() < tolerance
523    }
524
525    /// Construct the 3×3 rotation matrix corresponding to this unit quaternion.
526    ///
527    /// Returns the SO(3) matrix R such that R·v = `rotate_vector(v)` for all v ∈ ℝ³.
528    ///
529    /// # Lean Correspondence
530    ///
531    /// This is the Rust counterpart of `toRotationMatrix` in
532    /// `OrthogonalGroups.lean`. The Lean proofs establish:
533    /// - `toRotationMatrix_orthogonal_axiom`: R(q)ᵀ·R(q) = I
534    /// - `toRotationMatrix_det_axiom`: det(R(q)) = 1
535    /// - `toRotationMatrix_mul_axiom`: R(q₁·q₂) = R(q₁)·R(q₂)
536    /// - `toRotationMatrix_neg`: R(-q) = R(q)
537    #[must_use]
538    pub fn to_rotation_matrix(&self) -> [[f64; 3]; 3] {
539        let (w, x, y, z) = (self.w, self.x, self.y, self.z);
540        [
541            [
542                1.0 - 2.0 * (y * y + z * z),
543                2.0 * (x * y - w * z),
544                2.0 * (x * z + w * y),
545            ],
546            [
547                2.0 * (x * y + w * z),
548                1.0 - 2.0 * (x * x + z * z),
549                2.0 * (y * z - w * x),
550            ],
551            [
552                2.0 * (x * z - w * y),
553                2.0 * (y * z + w * x),
554                1.0 - 2.0 * (x * x + y * y),
555            ],
556        ]
557    }
558}
559
560// Quaternion multiplication: Hamilton product
561impl Mul for UnitQuaternion {
562    type Output = Self;
563
564    fn mul(self, rhs: Self) -> Self {
565        // (w₁ + x₁i + y₁j + z₁k)(w₂ + x₂i + y₂j + z₂k)
566        Self::new(
567            self.w * rhs.w - self.x * rhs.x - self.y * rhs.y - self.z * rhs.z,
568            self.w * rhs.x + self.x * rhs.w + self.y * rhs.z - self.z * rhs.y,
569            self.w * rhs.y - self.x * rhs.z + self.y * rhs.w + self.z * rhs.x,
570            self.w * rhs.z + self.x * rhs.y - self.y * rhs.x + self.z * rhs.w,
571        )
572    }
573}
574
575impl MulAssign for UnitQuaternion {
576    fn mul_assign(&mut self, rhs: Self) {
577        *self = *self * rhs;
578    }
579}
580
581#[cfg(test)]
582mod tests {
583    use super::*;
584    use std::f64::consts::PI;
585
586    #[test]
587    fn test_identity() {
588        let q = UnitQuaternion::identity();
589        assert_eq!(q.w, 1.0);
590        assert_eq!(q.x, 0.0);
591        assert_eq!(q.y, 0.0);
592        assert_eq!(q.z, 0.0);
593        assert!(q.verify_unit(1e-10));
594    }
595
596    #[test]
597    fn test_normalization() {
598        let q = UnitQuaternion::new(1.0, 2.0, 3.0, 4.0);
599        assert!((q.norm() - 1.0).abs() < 1e-10);
600        assert!(q.verify_unit(1e-10));
601    }
602
603    #[test]
604    fn test_axis_angle_roundtrip() {
605        let axis = [1.0, 2.0, 3.0];
606        let angle = PI / 3.0;
607
608        let q = UnitQuaternion::from_axis_angle(axis, angle);
609        let (axis_out, angle_out) = q.to_axis_angle();
610
611        // Normalize input axis for comparison
612        let axis_norm = (axis[0] * axis[0] + axis[1] * axis[1] + axis[2] * axis[2]).sqrt();
613        let axis_normalized = [
614            axis[0] / axis_norm,
615            axis[1] / axis_norm,
616            axis[2] / axis_norm,
617        ];
618
619        assert!((angle - angle_out).abs() < 1e-10);
620        for i in 0..3 {
621            assert!((axis_normalized[i] - axis_out[i]).abs() < 1e-10);
622        }
623
624        println!(
625            "✓ Axis-angle roundtrip: axis=[{:.3},{:.3},{:.3}], angle={:.3}",
626            axis_out[0], axis_out[1], axis_out[2], angle_out
627        );
628    }
629
630    #[test]
631    fn test_rotation_x() {
632        let q = UnitQuaternion::rotation_x(PI / 2.0);
633
634        // Should rotate (0,1,0) to (0,0,1)
635        let v = [0.0, 1.0, 0.0];
636        let rotated = q.rotate_vector(v);
637
638        assert!(rotated[0].abs() < 1e-10);
639        assert!(rotated[1].abs() < 1e-10);
640        assert!((rotated[2] - 1.0).abs() < 1e-10);
641
642        println!(
643            "✓ Rotation X by π/2: (0,1,0) → ({:.3},{:.3},{:.3})",
644            rotated[0], rotated[1], rotated[2]
645        );
646    }
647
648    #[test]
649    fn test_quaternion_multiplication() {
650        let q1 = UnitQuaternion::rotation_x(PI / 4.0);
651        let q2 = UnitQuaternion::rotation_y(PI / 4.0);
652
653        let q3 = q1 * q2;
654
655        // Should still be unit quaternion
656        assert!(q3.verify_unit(1e-10));
657
658        // Composition of rotations
659        let v = [1.0, 0.0, 0.0];
660        let rotated_composed = q3.rotate_vector(v);
661        let rotated_separate = q1.rotate_vector(q2.rotate_vector(v));
662
663        for i in 0..3 {
664            assert!((rotated_composed[i] - rotated_separate[i]).abs() < 1e-10);
665        }
666
667        println!("✓ Quaternion multiplication preserves unitarity and composition");
668    }
669
670    #[test]
671    fn test_inverse() {
672        let q = UnitQuaternion::rotation_z(PI / 3.0);
673        let q_inv = q.inverse();
674
675        let product = q * q_inv;
676
677        assert!((product.w - 1.0).abs() < 1e-10);
678        assert!(product.x.abs() < 1e-10);
679        assert!(product.y.abs() < 1e-10);
680        assert!(product.z.abs() < 1e-10);
681
682        println!("✓ Inverse: q · q⁻¹ = identity");
683    }
684
685    #[test]
686    fn test_distance_to_identity() {
687        let q_identity = UnitQuaternion::identity();
688        assert!(q_identity.distance_to_identity() < 1e-10);
689
690        let q_pi = UnitQuaternion::rotation_x(PI);
691        assert!((q_pi.distance_to_identity() - PI).abs() < 1e-10);
692
693        let q_half_pi = UnitQuaternion::rotation_y(PI / 2.0);
694        assert!((q_half_pi.distance_to_identity() - PI / 2.0).abs() < 1e-10);
695
696        println!(
697            "✓ Distance to identity: d(I)=0, d(π)={:.3}, d(π/2)={:.3}",
698            q_pi.distance_to_identity(),
699            q_half_pi.distance_to_identity()
700        );
701    }
702
703    #[test]
704    fn test_exp_log_roundtrip() {
705        let v = [0.5, 0.3, 0.2];
706
707        let q = UnitQuaternion::exp(v);
708        let v_out = q.log();
709
710        for i in 0..3 {
711            assert!((v[i] - v_out[i]).abs() < 1e-10);
712        }
713
714        println!(
715            "✓ Exp-log roundtrip: v=[{:.3},{:.3},{:.3}]",
716            v[0], v[1], v[2]
717        );
718    }
719
720    #[test]
721    fn test_slerp() {
722        let q1 = UnitQuaternion::rotation_x(0.0);
723        let q2 = UnitQuaternion::rotation_x(PI / 2.0);
724
725        let q_mid = q1.slerp(&q2, 0.5);
726
727        // Should be rotation by π/4
728        assert!((q_mid.distance_to_identity() - PI / 4.0).abs() < 1e-10);
729
730        // At t=0, should be q1
731        let q_0 = q1.slerp(&q2, 0.0);
732        assert!(q_0.distance_to(&q1) < 1e-10);
733
734        // At t=1, should be q2
735        let q_1 = q1.slerp(&q2, 1.0);
736        assert!(q_1.distance_to(&q2) < 1e-10);
737
738        println!("✓ SLERP: smooth interpolation on S³");
739    }
740
741    #[test]
742    fn test_near_pi_axis_angle_roundtrip() {
743        // Near θ = π, sin(θ/2) ≈ 1 and cos(θ/2) ≈ 0.
744        // Axis extraction divides by sin(θ/2), which is well-conditioned here.
745        // But w ≈ 0 means acos(w) sensitivity is low — this should work fine.
746        let axis = [0.0, 0.0, 1.0];
747        let angle = PI - 1e-8;
748
749        let q = UnitQuaternion::from_axis_angle(axis, angle);
750        assert!(q.verify_unit(1e-10));
751
752        let (axis_out, angle_out) = q.to_axis_angle();
753        assert!(
754            (angle - angle_out).abs() < 1e-6,
755            "Near-π angle roundtrip: got {}, expected {}",
756            angle_out,
757            angle
758        );
759        assert!(
760            (axis_out[2] - 1.0).abs() < 1e-6,
761            "Near-π axis roundtrip: got {:?}",
762            axis_out
763        );
764
765        // Also test the rotation acts correctly
766        let v = [1.0, 0.0, 0.0];
767        let rotated = q.rotate_vector(v);
768        // Rotation by π-ε around z should map (1,0,0) ≈ (-1, 0, 0)
769        assert!(
770            (rotated[0] + 1.0).abs() < 1e-6,
771            "Near-π rotation: got {:?}",
772            rotated
773        );
774    }
775
776    #[test]
777    fn test_exp_log_near_identity() {
778        // Very small algebra element — tests the near-identity branch
779        let v = [1e-12, 0.0, 0.0];
780        let q = UnitQuaternion::exp(v);
781        assert!(q.verify_unit(1e-14));
782        assert!(q.distance_to_identity() < 1e-10);
783
784        let log_q = q.log();
785        // Near identity, log should return near-zero
786        for &c in &log_q {
787            assert!(c.abs() < 1e-8);
788        }
789    }
790
791    #[test]
792    fn test_matrix_conversion() {
793        let q = UnitQuaternion::rotation_z(PI / 3.0);
794
795        let matrix = q.to_matrix();
796        let q_back = UnitQuaternion::from_matrix(matrix);
797
798        assert!((q.w - q_back.w).abs() < 1e-10);
799        assert!((q.x - q_back.x).abs() < 1e-10);
800        assert!((q.y - q_back.y).abs() < 1e-10);
801        assert!((q.z - q_back.z).abs() < 1e-10);
802
803        println!("✓ Matrix conversion roundtrip preserves quaternion");
804    }
805}