Skip to main content

lie_groups/
so3.rs

1//! Lie group SO(3) - 3D rotation group
2//!
3//! SO(3) is the group of 3×3 real orthogonal matrices with determinant 1.
4//! It represents rotations in 3-dimensional Euclidean space.
5//!
6//! # Mathematical Structure
7//!
8//! ```text
9//! SO(3) = { R ∈ ℝ³ˣ³ | R^T R = I, det(R) = 1 }
10//! ```
11//!
12//! # Lie Algebra
13//!
14//! The Lie algebra so(3) consists of 3×3 real antisymmetric matrices:
15//! ```text
16//! so(3) = { X ∈ ℝ³ˣ³ | X^T = -X }
17//! ```
18//!
19//! This is 3-dimensional, naturally isomorphic to ℝ³ via the cross product:
20//! ```text
21//! X ↔ v  where  X·w = v × w  for all w ∈ ℝ³
22//! ```
23//!
24//! # Relationship to SU(2)
25//!
26//! SO(3) ≅ SU(2)/ℤ₂ — SU(2) is the double cover of SO(3).
27//! Both have 3-dimensional Lie algebras: so(3) ≅ su(2) ≅ ℝ³
28//!
29//! # Sign Convention
30//!
31//! The so(3) bracket is the cross product: `[v, w] = v × w`, giving
32//! structure constants `fᵢⱼₖ = +εᵢⱼₖ`. The su(2) bracket (in the `{iσ/2}`
33//! basis) uses `[X, Y] = -(X × Y)`, giving `fᵢⱼₖ = -εᵢⱼₖ`.
34//!
35//! The bracket-preserving isomorphism φ: su(2) → so(3) maps `eₐ ↦ -Lₐ`,
36//! so that `φ([eₐ, eᵦ]) = [φ(eₐ), φ(eᵦ)]`. This sign arises from the
37//! standard `{iσ/2}` basis choice for su(2).
38
39use crate::traits::{AntiHermitianByConstruction, LieAlgebra, LieGroup, TracelessByConstruction};
40use nalgebra::Matrix3;
41use std::fmt;
42use std::ops::{Add, Mul, MulAssign, Neg, Sub};
43
44/// Lie algebra so(3) ≅ ℝ³
45///
46/// Elements of so(3) are 3×3 real antisymmetric matrices, which we represent
47/// as 3-vectors via the natural isomorphism with ℝ³.
48///
49/// # Isomorphism with ℝ³
50///
51/// An element v = (x, y, z) ∈ ℝ³ corresponds to the antisymmetric matrix:
52/// ```text
53/// [v]_× = [[0, -z, y], [z, 0, -x], [-y, x, 0]]
54/// ```
55///
56/// This satisfies: `[v]_×` · w = v × w (cross product)
57///
58/// # Basis Elements
59///
60/// The standard basis corresponds to angular momentum operators:
61/// ```text
62/// L_x = (1, 0, 0) ↔ [[0, 0, 0], [0, 0, -1], [0, 1, 0]]
63/// L_y = (0, 1, 0) ↔ [[0, 0, 1], [0, 0, 0], [-1, 0, 0]]
64/// L_z = (0, 0, 1) ↔ [[0, -1, 0], [1, 0, 0], [0, 0, 0]]
65/// ```
66#[derive(Clone, Copy, Debug, PartialEq)]
67pub struct So3Algebra(pub(crate) [f64; 3]);
68
69impl So3Algebra {
70    /// Create a new so(3) algebra element from components.
71    ///
72    /// The components `[x, y, z]` correspond to the antisymmetric matrix
73    /// `[[0, -z, y], [z, 0, -x], [-y, x, 0]]`.
74    #[must_use]
75    pub fn new(components: [f64; 3]) -> Self {
76        Self(components)
77    }
78
79    /// Returns the components as a fixed-size array reference.
80    #[must_use]
81    pub fn components(&self) -> &[f64; 3] {
82        &self.0
83    }
84}
85
86impl Add for So3Algebra {
87    type Output = Self;
88    fn add(self, rhs: Self) -> Self {
89        Self([
90            self.0[0] + rhs.0[0],
91            self.0[1] + rhs.0[1],
92            self.0[2] + rhs.0[2],
93        ])
94    }
95}
96
97impl Add<&So3Algebra> for So3Algebra {
98    type Output = So3Algebra;
99    fn add(self, rhs: &So3Algebra) -> So3Algebra {
100        self + *rhs
101    }
102}
103
104impl Add<So3Algebra> for &So3Algebra {
105    type Output = So3Algebra;
106    fn add(self, rhs: So3Algebra) -> So3Algebra {
107        *self + rhs
108    }
109}
110
111impl Add<&So3Algebra> for &So3Algebra {
112    type Output = So3Algebra;
113    fn add(self, rhs: &So3Algebra) -> So3Algebra {
114        *self + *rhs
115    }
116}
117
118impl Sub for So3Algebra {
119    type Output = Self;
120    fn sub(self, rhs: Self) -> Self {
121        Self([
122            self.0[0] - rhs.0[0],
123            self.0[1] - rhs.0[1],
124            self.0[2] - rhs.0[2],
125        ])
126    }
127}
128
129impl Neg for So3Algebra {
130    type Output = Self;
131    fn neg(self) -> Self {
132        Self([-self.0[0], -self.0[1], -self.0[2]])
133    }
134}
135
136impl Mul<f64> for So3Algebra {
137    type Output = Self;
138    fn mul(self, scalar: f64) -> Self {
139        Self([self.0[0] * scalar, self.0[1] * scalar, self.0[2] * scalar])
140    }
141}
142
143impl Mul<So3Algebra> for f64 {
144    type Output = So3Algebra;
145    fn mul(self, rhs: So3Algebra) -> So3Algebra {
146        rhs * self
147    }
148}
149
150impl LieAlgebra for So3Algebra {
151    const DIM: usize = 3;
152
153    #[inline]
154    fn zero() -> Self {
155        Self([0.0, 0.0, 0.0])
156    }
157
158    #[inline]
159    fn add(&self, other: &Self) -> Self {
160        Self([
161            self.0[0] + other.0[0],
162            self.0[1] + other.0[1],
163            self.0[2] + other.0[2],
164        ])
165    }
166
167    #[inline]
168    fn scale(&self, scalar: f64) -> Self {
169        Self([self.0[0] * scalar, self.0[1] * scalar, self.0[2] * scalar])
170    }
171
172    #[inline]
173    fn norm(&self) -> f64 {
174        (self.0[0].powi(2) + self.0[1].powi(2) + self.0[2].powi(2)).sqrt()
175    }
176
177    #[inline]
178    fn basis_element(i: usize) -> Self {
179        assert!(i < 3, "SO(3) algebra is 3-dimensional");
180        let mut v = [0.0; 3];
181        v[i] = 1.0;
182        Self(v)
183    }
184
185    #[inline]
186    fn from_components(components: &[f64]) -> Self {
187        assert_eq!(components.len(), 3, "so(3) has dimension 3");
188        Self([components[0], components[1], components[2]])
189    }
190
191    #[inline]
192    fn to_components(&self) -> Vec<f64> {
193        self.0.to_vec()
194    }
195
196    /// Lie bracket for so(3): [v, w] = v × w (cross product)
197    ///
198    /// For so(3) ≅ ℝ³, the Lie bracket is exactly the vector cross product.
199    ///
200    /// # Properties
201    ///
202    /// - Antisymmetric: [v, w] = -[w, v]
203    /// - Jacobi identity: [u, [v, w]] + [v, [w, u]] + [w, [u, v]] = 0
204    ///
205    /// # Examples
206    ///
207    /// ```
208    /// use lie_groups::so3::So3Algebra;
209    /// use lie_groups::traits::LieAlgebra;
210    ///
211    /// let lx = So3Algebra::basis_element(0);  // (1, 0, 0)
212    /// let ly = So3Algebra::basis_element(1);  // (0, 1, 0)
213    /// let bracket = lx.bracket(&ly);           // (1,0,0) × (0,1,0) = (0,0,1)
214    ///
215    /// // Should give L_z = (0, 0, 1)
216    /// assert!((bracket.components()[0]).abs() < 1e-10);
217    /// assert!((bracket.components()[1]).abs() < 1e-10);
218    /// assert!((bracket.components()[2] - 1.0).abs() < 1e-10);
219    /// ```
220    #[inline]
221    fn bracket(&self, other: &Self) -> Self {
222        // Cross product: v × w
223        let v = self.0;
224        let w = other.0;
225
226        Self([
227            v[1] * w[2] - v[2] * w[1], // x component
228            v[2] * w[0] - v[0] * w[2], // y component
229            v[0] * w[1] - v[1] * w[0], // z component
230        ])
231    }
232
233    #[inline]
234    fn inner(&self, other: &Self) -> f64 {
235        self.0[0] * other.0[0] + self.0[1] * other.0[1] + self.0[2] * other.0[2]
236    }
237}
238
239/// SO(3) group element - 3×3 real orthogonal matrix with determinant 1
240///
241/// Represents a rotation in 3D space.
242///
243/// # Representation
244///
245/// We use `Matrix3<f64>` from nalgebra to represent the 3×3 rotation matrix.
246///
247/// # Constraints
248///
249/// - Orthogonality: R^T R = I
250/// - Determinant: det(R) = 1
251///
252/// # Examples
253///
254/// ```
255/// use lie_groups::so3::SO3;
256/// use lie_groups::traits::LieGroup;
257///
258/// // Rotation around Z-axis by π/2
259/// let rot = SO3::rotation_z(std::f64::consts::FRAC_PI_2);
260///
261/// // Verify it's orthogonal
262/// assert!(rot.verify_orthogonality(1e-10));
263/// ```
264#[derive(Clone, Debug, PartialEq)]
265pub struct SO3 {
266    /// 3×3 real orthogonal matrix
267    pub(crate) matrix: Matrix3<f64>,
268}
269
270impl SO3 {
271    /// Access the underlying 3×3 orthogonal matrix
272    #[must_use]
273    pub fn matrix(&self) -> &Matrix3<f64> {
274        &self.matrix
275    }
276
277    /// Identity element (no rotation)
278    #[must_use]
279    pub fn identity() -> Self {
280        Self {
281            matrix: Matrix3::identity(),
282        }
283    }
284
285    /// Rotation around X-axis by angle θ (in radians)
286    ///
287    /// ```text
288    /// R_x(θ) = [[1, 0, 0], [0, cos(θ), -sin(θ)], [0, sin(θ), cos(θ)]]
289    /// ```
290    #[must_use]
291    pub fn rotation_x(angle: f64) -> Self {
292        let (s, c) = angle.sin_cos();
293
294        Self {
295            matrix: Matrix3::new(1.0, 0.0, 0.0, 0.0, c, -s, 0.0, s, c),
296        }
297    }
298
299    /// Rotation around Y-axis by angle θ (in radians)
300    ///
301    /// ```text
302    /// R_y(θ) = [[cos(θ), 0, sin(θ)], [0, 1, 0], [-sin(θ), 0, cos(θ)]]
303    /// ```
304    #[must_use]
305    pub fn rotation_y(angle: f64) -> Self {
306        let (s, c) = angle.sin_cos();
307
308        Self {
309            matrix: Matrix3::new(c, 0.0, s, 0.0, 1.0, 0.0, -s, 0.0, c),
310        }
311    }
312
313    /// Rotation around Z-axis by angle θ (in radians)
314    ///
315    /// ```text
316    /// R_z(θ) = [[cos(θ), -sin(θ), 0], [sin(θ), cos(θ), 0], [0, 0, 1]]
317    /// ```
318    #[must_use]
319    pub fn rotation_z(angle: f64) -> Self {
320        let (s, c) = angle.sin_cos();
321
322        Self {
323            matrix: Matrix3::new(c, -s, 0.0, s, c, 0.0, 0.0, 0.0, 1.0),
324        }
325    }
326
327    /// Rotation around arbitrary axis by angle
328    ///
329    /// Uses Rodrigues' rotation formula:
330    /// ```text
331    /// R(θ, n̂) = I + sin(θ)[n̂]_× + (1-cos(θ))[n̂]_ײ
332    /// ```
333    ///
334    /// where [n̂]_× is the skew-symmetric matrix for axis n̂.
335    #[must_use]
336    pub fn rotation(axis: [f64; 3], angle: f64) -> Self {
337        let norm = (axis[0].powi(2) + axis[1].powi(2) + axis[2].powi(2)).sqrt();
338        if norm < 1e-10 {
339            return Self::identity();
340        }
341
342        // Normalize axis
343        let n = [axis[0] / norm, axis[1] / norm, axis[2] / norm];
344
345        // Rodrigues formula
346        let (s, c) = angle.sin_cos();
347        let t = 1.0 - c;
348
349        let matrix = Matrix3::new(
350            t * n[0] * n[0] + c,
351            t * n[0] * n[1] - s * n[2],
352            t * n[0] * n[2] + s * n[1],
353            t * n[0] * n[1] + s * n[2],
354            t * n[1] * n[1] + c,
355            t * n[1] * n[2] - s * n[0],
356            t * n[0] * n[2] - s * n[1],
357            t * n[1] * n[2] + s * n[0],
358            t * n[2] * n[2] + c,
359        );
360
361        Self { matrix }
362    }
363
364    /// Trace of the rotation matrix: Tr(R) = 1 + 2cos(θ)
365    #[must_use]
366    pub fn trace(&self) -> f64 {
367        self.matrix.trace()
368    }
369
370    /// Convert to 3×3 array format
371    #[must_use]
372    pub fn to_matrix(&self) -> [[f64; 3]; 3] {
373        let m = &self.matrix;
374        [
375            [m[(0, 0)], m[(0, 1)], m[(0, 2)]],
376            [m[(1, 0)], m[(1, 1)], m[(1, 2)]],
377            [m[(2, 0)], m[(2, 1)], m[(2, 2)]],
378        ]
379    }
380
381    /// Create from 3×3 array format
382    ///
383    /// # Panics
384    ///
385    /// Does not verify orthogonality. Use `verify_orthogonality()` after construction.
386    #[must_use]
387    pub fn from_matrix(arr: [[f64; 3]; 3]) -> Self {
388        Self {
389            matrix: Matrix3::new(
390                arr[0][0], arr[0][1], arr[0][2], arr[1][0], arr[1][1], arr[1][2], arr[2][0],
391                arr[2][1], arr[2][2],
392            ),
393        }
394    }
395
396    /// Verify orthogonality: R^T R = I
397    #[must_use]
398    pub fn verify_orthogonality(&self, tolerance: f64) -> bool {
399        let product = self.matrix.transpose() * self.matrix;
400        let identity = Matrix3::identity();
401
402        (product - identity).norm() < tolerance
403    }
404
405    /// Matrix inverse (equals transpose for orthogonal matrices)
406    #[must_use]
407    pub fn inverse(&self) -> Self {
408        Self {
409            matrix: self.matrix.transpose(),
410        }
411    }
412
413    /// Distance from identity (rotation angle)
414    ///
415    /// For a rotation matrix R, the angle θ satisfies:
416    /// ```text
417    /// trace(R) = 1 + 2cos(θ)
418    /// ```
419    #[must_use]
420    pub fn distance_to_identity(&self) -> f64 {
421        let trace = self.matrix.trace();
422        let cos_theta = (trace - 1.0) / 2.0;
423
424        // Clamp to [-1, 1] to handle numerical errors
425        let cos_theta = cos_theta.clamp(-1.0, 1.0);
426
427        cos_theta.acos()
428    }
429
430    /// Interpolate between two SO(3) elements with proper orthogonalization
431    ///
432    /// Uses linear interpolation of matrix elements followed by Gram-Schmidt
433    /// orthogonalization to ensure the result stays on SO(3).
434    ///
435    /// # Arguments
436    /// * `other` - The target rotation
437    /// * `t` - Interpolation parameter in [0, 1]
438    ///
439    /// # Returns
440    /// An SO(3) element: self at t=0, other at t=1
441    ///
442    /// # Note
443    /// This is NOT geodesic interpolation (SLERP). For true geodesic paths,
444    /// convert to quaternions and use quaternion SLERP. However, this method
445    /// guarantees the result is always a valid rotation matrix.
446    #[must_use]
447    pub fn interpolate(&self, other: &Self, t: f64) -> Self {
448        if t <= 0.0 {
449            return self.clone();
450        }
451        if t >= 1.0 {
452            return other.clone();
453        }
454
455        // Linear interpolation of matrix elements
456        let interpolated = self.matrix * (1.0 - t) + other.matrix * t;
457
458        // Gram-Schmidt orthogonalization to project back onto SO(3)
459        Self::gram_schmidt_orthogonalize(interpolated)
460    }
461
462    /// Orthogonalize a matrix using Gram-Schmidt process
463    ///
464    /// Takes a near-orthogonal matrix and projects it back onto SO(3).
465    /// This is essential for numerical stability when matrices drift
466    /// from orthogonality due to floating-point accumulation.
467    ///
468    /// # Algorithm
469    /// 1. Normalize first column
470    /// 2. Orthogonalize and normalize second column
471    /// 3. Compute third column as cross product (ensures determinant = 1)
472    #[must_use]
473    pub fn gram_schmidt_orthogonalize(matrix: Matrix3<f64>) -> Self {
474        use nalgebra::Vector3;
475
476        // Extract columns
477        let c0 = Vector3::new(matrix[(0, 0)], matrix[(1, 0)], matrix[(2, 0)]);
478        let c1 = Vector3::new(matrix[(0, 1)], matrix[(1, 1)], matrix[(2, 1)]);
479
480        // Gram-Schmidt orthogonalization
481        let e0 = c0.normalize();
482        let e1 = (c1 - e0 * e0.dot(&c1)).normalize();
483        // Third column is cross product to ensure det = +1 (proper rotation)
484        let e2 = e0.cross(&e1);
485
486        let orthogonal = Matrix3::new(
487            e0[0], e1[0], e2[0], e0[1], e1[1], e2[1], e0[2], e1[2], e2[2],
488        );
489
490        Self { matrix: orthogonal }
491    }
492
493    /// Re-normalize a rotation matrix that may have drifted
494    ///
495    /// Use this periodically after many matrix multiplications to
496    /// prevent numerical drift from accumulating.
497    #[must_use]
498    pub fn renormalize(&self) -> Self {
499        Self::gram_schmidt_orthogonalize(self.matrix)
500    }
501
502    /// Geodesic distance between two SO(3) elements
503    ///
504    /// d(R₁, R₂) = ||log(R₁ᵀ R₂)||
505    ///
506    /// This is the rotation angle of the relative rotation R₁ᵀ R₂.
507    #[must_use]
508    pub fn geodesic_distance(&self, other: &Self) -> f64 {
509        // Compute relative rotation: R_rel = R₁ᵀ R₂
510        let relative = Self {
511            matrix: self.matrix.transpose() * other.matrix,
512        };
513        relative.distance_to_identity()
514    }
515}
516
517impl approx::AbsDiffEq for So3Algebra {
518    type Epsilon = f64;
519
520    fn default_epsilon() -> Self::Epsilon {
521        1e-10
522    }
523
524    fn abs_diff_eq(&self, other: &Self, epsilon: Self::Epsilon) -> bool {
525        self.0
526            .iter()
527            .zip(other.0.iter())
528            .all(|(a, b)| (a - b).abs() < epsilon)
529    }
530}
531
532impl approx::RelativeEq for So3Algebra {
533    fn default_max_relative() -> Self::Epsilon {
534        1e-10
535    }
536
537    fn relative_eq(
538        &self,
539        other: &Self,
540        epsilon: Self::Epsilon,
541        max_relative: Self::Epsilon,
542    ) -> bool {
543        self.0
544            .iter()
545            .zip(other.0.iter())
546            .all(|(a, b)| approx::RelativeEq::relative_eq(a, b, epsilon, max_relative))
547    }
548}
549
550impl fmt::Display for So3Algebra {
551    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
552        write!(
553            f,
554            "so(3)[{:.4}, {:.4}, {:.4}]",
555            self.0[0], self.0[1], self.0[2]
556        )
557    }
558}
559
560impl fmt::Display for SO3 {
561    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
562        let dist = self.distance_to_identity();
563        write!(f, "SO(3)(θ={:.4})", dist)
564    }
565}
566
567/// Group multiplication: R₁ · R₂
568impl Mul<&SO3> for &SO3 {
569    type Output = SO3;
570    fn mul(self, rhs: &SO3) -> SO3 {
571        SO3 {
572            matrix: self.matrix * rhs.matrix,
573        }
574    }
575}
576
577impl Mul<&SO3> for SO3 {
578    type Output = SO3;
579    fn mul(self, rhs: &SO3) -> SO3 {
580        &self * rhs
581    }
582}
583
584impl MulAssign<&SO3> for SO3 {
585    fn mul_assign(&mut self, rhs: &SO3) {
586        self.matrix *= rhs.matrix;
587    }
588}
589
590impl LieGroup for SO3 {
591    const MATRIX_DIM: usize = 3;
592
593    type Algebra = So3Algebra;
594
595    fn identity() -> Self {
596        Self::identity()
597    }
598
599    fn compose(&self, other: &Self) -> Self {
600        Self {
601            matrix: self.matrix * other.matrix,
602        }
603    }
604
605    fn inverse(&self) -> Self {
606        Self::inverse(self)
607    }
608
609    fn conjugate_transpose(&self) -> Self {
610        // For orthogonal matrices, conjugate transpose = transpose = inverse
611        self.inverse()
612    }
613
614    fn adjoint_action(&self, algebra_element: &So3Algebra) -> So3Algebra {
615        // For SO(3): Ad_R(v) = R v (matrix-vector multiplication)
616        //
617        // This corresponds to rotating the axis v by rotation R
618        let v = algebra_element.0;
619        let rotated = self.matrix * nalgebra::Vector3::new(v[0], v[1], v[2]);
620
621        So3Algebra([rotated[0], rotated[1], rotated[2]])
622    }
623
624    fn distance_to_identity(&self) -> f64 {
625        Self::distance_to_identity(self)
626    }
627
628    fn exp(tangent: &So3Algebra) -> Self {
629        // Exponential map using Rodrigues formula
630        let angle = tangent.norm();
631
632        if angle < 1e-10 {
633            return Self::identity();
634        }
635
636        let axis = [
637            tangent.0[0] / angle,
638            tangent.0[1] / angle,
639            tangent.0[2] / angle,
640        ];
641
642        Self::rotation(axis, angle)
643    }
644
645    fn log(&self) -> crate::error::LogResult<So3Algebra> {
646        use crate::error::LogError;
647
648        // For SO(3), extract rotation angle and axis from the rotation matrix
649        //
650        // The angle θ satisfies:
651        // Tr(R) = 1 + 2·cos(θ)
652        //
653        // The axis n̂ can be extracted from the skew-symmetric part:
654        // R - R^T = 2·sin(θ)·[n̂]×
655        // where [n̂]× is the cross-product matrix
656
657        let trace = self.matrix.trace();
658
659        // Extract angle: θ = arccos((Tr(R) - 1) / 2)
660        let cos_theta = ((trace - 1.0) / 2.0).clamp(-1.0, 1.0);
661        let theta = cos_theta.acos();
662
663        const SMALL_ANGLE_THRESHOLD: f64 = 1e-10;
664        const SINGULARITY_THRESHOLD: f64 = 1e-6;
665
666        if theta.abs() < SMALL_ANGLE_THRESHOLD {
667            // Near identity: log(I) = 0
668            return Ok(So3Algebra::zero());
669        }
670
671        if (theta - std::f64::consts::PI).abs() < SINGULARITY_THRESHOLD {
672            // Near singularity at θ = π (rotation by 180°)
673            // At exactly θ = π, the axis direction is ambiguous
674            return Err(LogError::Singularity {
675                reason: format!(
676                    "SO(3) element at rotation angle θ ≈ π (θ = {:.6}), axis is ambiguous",
677                    theta
678                ),
679            });
680        }
681
682        // Extract rotation axis from skew-symmetric part
683        // R - R^T = 2·sin(θ)·[n̂]×
684        //
685        // [n̂]× = [[ 0,  -nz,  ny],
686        //         [ nz,  0,  -nx],
687        //         [-ny,  nx,  0 ]]
688        //
689        // So:
690        // R[2,1] - R[1,2] = 2·sin(θ)·nx
691        // R[0,2] - R[2,0] = 2·sin(θ)·ny
692        // R[1,0] - R[0,1] = 2·sin(θ)·nz
693
694        let sin_theta = theta.sin();
695
696        if sin_theta.abs() < 1e-10 {
697            // This shouldn't happen given our checks above, but guard against it
698            return Err(LogError::NumericalInstability {
699                reason: "sin(θ) too small for reliable axis extraction".to_string(),
700            });
701        }
702
703        let two_sin_theta = 2.0 * sin_theta;
704        let nx = (self.matrix[(2, 1)] - self.matrix[(1, 2)]) / two_sin_theta;
705        let ny = (self.matrix[(0, 2)] - self.matrix[(2, 0)]) / two_sin_theta;
706        let nz = (self.matrix[(1, 0)] - self.matrix[(0, 1)]) / two_sin_theta;
707
708        // The logarithm is log(R) = θ·(nx, ny, nz) ∈ so(3)
709        Ok(So3Algebra([theta * nx, theta * ny, theta * nz]))
710    }
711}
712
713// ============================================================================
714// Mathematical Property Implementations
715// ============================================================================
716
717use crate::traits::{Compact, SemiSimple, Simple};
718
719/// SO(3) is compact
720///
721/// The rotation group is diffeomorphic to ℝP³ (real projective 3-space).
722/// All rotations are bounded: ||R|| = 1.
723impl Compact for SO3 {}
724
725/// SO(3) is simple
726///
727/// It has no non-trivial normal subgroups (for dimension > 2).
728impl Simple for SO3 {}
729
730/// SO(3) is semi-simple
731impl SemiSimple for SO3 {}
732
733// ============================================================================
734// Algebra Marker Traits
735// ============================================================================
736
737/// so(3) algebra elements are traceless by construction.
738///
739/// The representation `So3Algebra::new([f64; 3])` stores coefficients for
740/// 3×3 antisymmetric matrices. All antisymmetric matrices are traceless.
741impl TracelessByConstruction for So3Algebra {}
742
743/// so(3) algebra elements are anti-Hermitian by construction.
744///
745/// Real antisymmetric matrices satisfy A^T = -A, which over ℝ is
746/// equivalent to A† = -A (anti-Hermitian).
747impl AntiHermitianByConstruction for So3Algebra {}
748
749#[cfg(test)]
750mod tests {
751    use super::*;
752    use approx::assert_relative_eq;
753
754    #[test]
755    fn test_identity() {
756        let id = SO3::identity();
757        assert!(id.verify_orthogonality(1e-10));
758        assert_relative_eq!(id.distance_to_identity(), 0.0, epsilon = 1e-10);
759    }
760
761    #[test]
762    fn test_rotations_orthogonal() {
763        let rx = SO3::rotation_x(0.5);
764        let ry = SO3::rotation_y(1.2);
765        let rz = SO3::rotation_z(2.1);
766
767        assert!(rx.verify_orthogonality(1e-10));
768        assert!(ry.verify_orthogonality(1e-10));
769        assert!(rz.verify_orthogonality(1e-10));
770    }
771
772    #[test]
773    fn test_inverse() {
774        let r = SO3::rotation_x(0.7);
775        let r_inv = r.inverse();
776        let product = r.compose(&r_inv);
777
778        assert_relative_eq!(product.distance_to_identity(), 0.0, epsilon = 1e-10);
779    }
780
781    #[test]
782    fn test_bracket_antisymmetry() {
783        use crate::traits::LieAlgebra;
784
785        let v = So3Algebra([0.7, -0.3, 1.2]);
786        let w = So3Algebra([-0.5, 0.9, 0.4]);
787
788        let vw = v.bracket(&w);
789        let wv = w.bracket(&v);
790
791        // [v, w] = -[w, v]
792        for i in 0..3 {
793            assert_relative_eq!(vw.0[i], -wv.0[i], epsilon = 1e-14);
794        }
795    }
796
797    #[test]
798    fn test_bracket_bilinearity() {
799        use crate::traits::LieAlgebra;
800
801        let x = So3Algebra([1.0, 0.0, 0.0]);
802        let y = So3Algebra([0.0, 1.0, 0.0]);
803        let z = So3Algebra([0.0, 0.0, 1.0]);
804        let alpha = 2.5;
805
806        // Left linearity: [αX + Y, Z] = α[X, Z] + [Y, Z]
807        let lhs = x.scale(alpha).add(&y).bracket(&z);
808        let rhs = x.bracket(&z).scale(alpha).add(&y.bracket(&z));
809        for i in 0..3 {
810            assert_relative_eq!(lhs.0[i], rhs.0[i], epsilon = 1e-14);
811        }
812
813        // Right linearity: [Z, αX + Y] = α[Z, X] + [Z, Y]
814        let lhs = z.bracket(&x.scale(alpha).add(&y));
815        let rhs = z.bracket(&x).scale(alpha).add(&z.bracket(&y));
816        for i in 0..3 {
817            assert_relative_eq!(lhs.0[i], rhs.0[i], epsilon = 1e-14);
818        }
819    }
820
821    #[test]
822    fn test_algebra_bracket() {
823        use crate::traits::LieAlgebra;
824
825        let lx = So3Algebra::basis_element(0);
826        let ly = So3Algebra::basis_element(1);
827        let bracket = lx.bracket(&ly);
828
829        // [L_x, L_y] = L_z
830        assert_relative_eq!(bracket.0[0], 0.0, epsilon = 1e-10);
831        assert_relative_eq!(bracket.0[1], 0.0, epsilon = 1e-10);
832        assert_relative_eq!(bracket.0[2], 1.0, epsilon = 1e-10);
833    }
834
835    #[test]
836    fn test_adjoint_action() {
837        use crate::traits::LieGroup;
838
839        // Rotate L_x by 90° around Z
840        let rz = SO3::rotation_z(std::f64::consts::FRAC_PI_2);
841        let lx = So3Algebra([1.0, 0.0, 0.0]);
842        let rotated = rz.adjoint_action(&lx);
843
844        // Should give L_y = (0, 1, 0)
845        assert_relative_eq!(rotated.0[0], 0.0, epsilon = 1e-10);
846        assert_relative_eq!(rotated.0[1], 1.0, epsilon = 1e-10);
847        assert_relative_eq!(rotated.0[2], 0.0, epsilon = 1e-10);
848    }
849
850    #[test]
851    fn test_interpolate_orthogonality() {
852        // Interpolation should always produce orthogonal matrices
853        let r1 = SO3::rotation_x(0.3);
854        let r2 = SO3::rotation_y(1.2);
855
856        for t in [0.0, 0.25, 0.5, 0.75, 1.0] {
857            let interp = r1.interpolate(&r2, t);
858            assert!(
859                interp.verify_orthogonality(1e-10),
860                "Interpolated matrix not orthogonal at t={}",
861                t
862            );
863        }
864    }
865
866    #[test]
867    fn test_interpolate_endpoints() {
868        let r1 = SO3::rotation_z(0.5);
869        let r2 = SO3::rotation_z(1.5);
870
871        let at_0 = r1.interpolate(&r2, 0.0);
872        let at_1 = r1.interpolate(&r2, 1.0);
873
874        assert_relative_eq!(
875            at_0.distance_to_identity(),
876            r1.distance_to_identity(),
877            epsilon = 1e-10
878        );
879        assert_relative_eq!(
880            at_1.distance_to_identity(),
881            r2.distance_to_identity(),
882            epsilon = 1e-10
883        );
884    }
885
886    #[test]
887    fn test_gram_schmidt_orthogonalize() {
888        // Create a slightly perturbed (non-orthogonal) matrix
889        let perturbed = Matrix3::new(1.001, 0.01, 0.0, 0.0, 0.999, 0.01, 0.0, 0.0, 1.002);
890
891        let fixed = SO3::gram_schmidt_orthogonalize(perturbed);
892        assert!(fixed.verify_orthogonality(1e-10));
893
894        // Check determinant is +1 (proper rotation)
895        let det = fixed.matrix.determinant();
896        assert_relative_eq!(det, 1.0, epsilon = 1e-10);
897    }
898
899    #[test]
900    fn test_renormalize() {
901        // Simulate drift from many multiplications
902        let r = SO3::rotation_z(0.001);
903        let mut accumulated = SO3::identity();
904
905        for _ in 0..1000 {
906            accumulated = accumulated.compose(&r);
907        }
908
909        // May have drifted slightly
910        let renormalized = accumulated.renormalize();
911        assert!(renormalized.verify_orthogonality(1e-12));
912    }
913
914    #[test]
915    fn test_geodesic_distance_symmetric() {
916        let r1 = SO3::rotation_x(0.5);
917        let r2 = SO3::rotation_y(1.0);
918
919        let d12 = r1.geodesic_distance(&r2);
920        let d21 = r2.geodesic_distance(&r1);
921
922        assert_relative_eq!(d12, d21, epsilon = 1e-10);
923    }
924
925    #[test]
926    fn test_exp_log_roundtrip() {
927        use crate::traits::{LieAlgebra, LieGroup};
928
929        // Small angles (well within convergence)
930        for &angle in &[0.1, 0.5, 1.0, 2.0] {
931            let v = So3Algebra([angle * 0.6, angle * 0.8, 0.0]);
932            let r = SO3::exp(&v);
933            assert!(r.verify_orthogonality(1e-10));
934
935            let v_back = SO3::log(&r).expect("log should succeed for small angles");
936            let diff = v.add(&v_back.scale(-1.0));
937            assert!(
938                diff.norm() < 1e-9,
939                "exp/log roundtrip failed at angle {}: error = {:.2e}",
940                angle,
941                diff.norm()
942            );
943        }
944    }
945
946    #[test]
947    fn test_log_exp_roundtrip() {
948        use crate::traits::LieGroup;
949
950        // Start from group elements, take log, then exp back
951        let rotations = [
952            SO3::rotation_x(0.7),
953            SO3::rotation_y(1.3),
954            SO3::rotation_z(0.4),
955            SO3::rotation_x(0.3).compose(&SO3::rotation_y(0.5)),
956        ];
957
958        for (i, r) in rotations.iter().enumerate() {
959            let v = SO3::log(r).expect("log should succeed");
960            let r_back = SO3::exp(&v);
961            let dist = r.geodesic_distance(&r_back);
962            assert!(
963                dist < 1e-9,
964                "log/exp roundtrip failed for rotation {}: distance = {:.2e}",
965                i,
966                dist
967            );
968        }
969    }
970
971    #[test]
972    fn test_jacobi_identity() {
973        use crate::traits::LieAlgebra;
974
975        let x = So3Algebra([1.0, 0.0, 0.0]);
976        let y = So3Algebra([0.0, 1.0, 0.0]);
977        let z = So3Algebra([0.0, 0.0, 1.0]);
978
979        // [X, [Y, Z]] + [Y, [Z, X]] + [Z, [X, Y]] = 0
980        let t1 = x.bracket(&y.bracket(&z));
981        let t2 = y.bracket(&z.bracket(&x));
982        let t3 = z.bracket(&x.bracket(&y));
983        let sum = t1.add(&t2).add(&t3);
984
985        assert!(
986            sum.norm() < 1e-14,
987            "Jacobi identity violated: ||sum|| = {:.2e}",
988            sum.norm()
989        );
990    }
991
992    #[test]
993    fn test_geodesic_distance_to_self() {
994        let r = SO3::rotation([1.0, 2.0, 3.0], 0.8);
995        let d = r.geodesic_distance(&r);
996        assert_relative_eq!(d, 0.0, epsilon = 1e-10);
997    }
998
999    #[test]
1000    fn test_geodesic_distance_to_identity() {
1001        let angle = 0.7;
1002        let r = SO3::rotation_z(angle);
1003
1004        let d = SO3::identity().geodesic_distance(&r);
1005        assert_relative_eq!(d, angle, epsilon = 1e-10);
1006    }
1007
1008    // ========================================================================
1009    // Property-Based Tests for Group Axioms
1010    // ========================================================================
1011    //
1012    // These tests use proptest to verify that SO(3) satisfies the
1013    // mathematical axioms of a Lie group for randomly generated elements.
1014    //
1015    // This is a form of **specification-based testing**: the group axioms
1016    // are the specification, and we verify they hold for all inputs.
1017    //
1018    // Run with: cargo test --features nightly
1019
1020    #[cfg(feature = "proptest")]
1021    use proptest::prelude::*;
1022
1023    /// Strategy for generating arbitrary SO(3) elements.
1024    ///
1025    /// We generate SO(3) elements by composing three Euler rotations:
1026    /// `R = R_z(α) · R_y(β) · R_x(γ)`
1027    ///
1028    /// This gives good coverage of SO(3) ≅ ℝP³.
1029    #[cfg(feature = "proptest")]
1030    fn arb_so3() -> impl Strategy<Value = SO3> {
1031        use std::f64::consts::TAU;
1032
1033        // Generate three Euler angles
1034        let alpha = 0.0..TAU;
1035        let beta = 0.0..TAU;
1036        let gamma = 0.0..TAU;
1037
1038        (alpha, beta, gamma).prop_map(|(a, b, c)| {
1039            SO3::rotation_z(a)
1040                .compose(&SO3::rotation_y(b))
1041                .compose(&SO3::rotation_x(c))
1042        })
1043    }
1044
1045    /// Strategy for generating arbitrary `So3Algebra` elements.
1046    ///
1047    /// We generate algebra elements by picking random coefficients in [-π, π]
1048    /// for each of the three basis directions (`L_x`, `L_y`, `L_z`).
1049    #[cfg(feature = "proptest")]
1050    fn arb_so3_algebra() -> impl Strategy<Value = So3Algebra> {
1051        use std::f64::consts::PI;
1052
1053        ((-PI..PI), (-PI..PI), (-PI..PI)).prop_map(|(a, b, c)| So3Algebra([a, b, c]))
1054    }
1055
1056    #[cfg(feature = "proptest")]
1057    proptest! {
1058        /// **Group Axiom 1: Identity Element**
1059        ///
1060        /// For all R ∈ SO(3):
1061        /// - I · R = R (left identity)
1062        /// - R · I = R (right identity)
1063        ///
1064        /// where I = identity rotation
1065        ///
1066        /// Note: We use tolerance 1e-7 to account for floating-point
1067        /// rounding errors in matrix operations.
1068        #[test]
1069        fn prop_identity_axiom(r in arb_so3()) {
1070            let e = SO3::identity();
1071
1072            // Left identity: I · R = R
1073            let left = e.compose(&r);
1074            prop_assert!(
1075                left.distance(&r) < 1e-7,
1076                "Left identity failed: I·R != R, distance = {}",
1077                left.distance(&r)
1078            );
1079
1080            // Right identity: R · I = R
1081            let right = r.compose(&e);
1082            prop_assert!(
1083                right.distance(&r) < 1e-7,
1084                "Right identity failed: R·I != R, distance = {}",
1085                right.distance(&r)
1086            );
1087        }
1088
1089        /// **Group Axiom 2: Inverse Element**
1090        ///
1091        /// For all R ∈ SO(3):
1092        /// - R · R⁻¹ = I (right inverse)
1093        /// - R⁻¹ · R = I (left inverse)
1094        ///
1095        /// where R⁻¹ = inverse of R (equals R^T for orthogonal matrices)
1096        ///
1097        /// Note: We use tolerance 1e-7 to account for floating-point
1098        /// rounding errors in matrix operations.
1099        #[test]
1100        fn prop_inverse_axiom(r in arb_so3()) {
1101            let r_inv = r.inverse();
1102
1103            // Right inverse: R · R⁻¹ = I
1104            let right_product = r.compose(&r_inv);
1105            prop_assert!(
1106                right_product.is_near_identity(1e-7),
1107                "Right inverse failed: R·R⁻¹ != I, distance = {}",
1108                right_product.distance_to_identity()
1109            );
1110
1111            // Left inverse: R⁻¹ · R = I
1112            let left_product = r_inv.compose(&r);
1113            prop_assert!(
1114                left_product.is_near_identity(1e-7),
1115                "Left inverse failed: R⁻¹·R != I, distance = {}",
1116                left_product.distance_to_identity()
1117            );
1118        }
1119
1120        /// **Group Axiom 3: Associativity**
1121        ///
1122        /// For all R₁, R₂, R₃ ∈ SO(3):
1123        /// - (R₁ · R₂) · R₃ = R₁ · (R₂ · R₃)
1124        ///
1125        /// Group multiplication is associative.
1126        ///
1127        /// Note: We use tolerance 1e-7 to account for floating-point
1128        /// rounding errors in matrix operations.
1129        #[test]
1130        fn prop_associativity(r1 in arb_so3(), r2 in arb_so3(), r3 in arb_so3()) {
1131            // Left association: (R₁ · R₂) · R₃
1132            let left_assoc = r1.compose(&r2).compose(&r3);
1133
1134            // Right association: R₁ · (R₂ · R₃)
1135            let right_assoc = r1.compose(&r2.compose(&r3));
1136
1137            prop_assert!(
1138                left_assoc.distance(&right_assoc) < 1e-7,
1139                "Associativity failed: (R₁·R₂)·R₃ != R₁·(R₂·R₃), distance = {}",
1140                left_assoc.distance(&right_assoc)
1141            );
1142        }
1143
1144        /// **Lie Group Property: Inverse is Smooth**
1145        ///
1146        /// For SO(3), the inverse operation is smooth (continuously differentiable).
1147        /// We verify this by checking that nearby elements have nearby inverses.
1148        #[test]
1149        fn prop_inverse_continuity(r in arb_so3()) {
1150            // Create a small perturbation
1151            let epsilon = 0.01;
1152            let perturbation = SO3::rotation_x(epsilon);
1153            let r_perturbed = r.compose(&perturbation);
1154
1155            // Check that inverses are close
1156            let inv_distance = r.inverse().distance(&r_perturbed.inverse());
1157
1158            prop_assert!(
1159                inv_distance < 0.1,
1160                "Inverse not continuous: small perturbation caused large inverse change, distance = {}",
1161                inv_distance
1162            );
1163        }
1164
1165        /// **Orthogonality Preservation**
1166        ///
1167        /// All SO(3) operations should preserve orthogonality.
1168        /// This is not strictly a group axiom, but it's essential for SO(3).
1169        #[test]
1170        fn prop_orthogonality_preserved(r1 in arb_so3(), r2 in arb_so3()) {
1171            // Composition preserves orthogonality
1172            let product = r1.compose(&r2);
1173            prop_assert!(
1174                product.verify_orthogonality(1e-10),
1175                "Composition violated orthogonality"
1176            );
1177
1178            // Inverse preserves orthogonality
1179            let inv = r1.inverse();
1180            prop_assert!(
1181                inv.verify_orthogonality(1e-10),
1182                "Inverse violated orthogonality"
1183            );
1184        }
1185
1186        /// **Adjoint Representation: Group Homomorphism**
1187        ///
1188        /// The adjoint representation Ad: G → Aut(𝔤) is a group homomorphism:
1189        /// - Ad_{R₁∘R₂}(v) = Ad_{R₁}(Ad_{R₂}(v))
1190        ///
1191        /// This is a fundamental property that must hold for the adjoint action
1192        /// to be a valid representation of the group.
1193        #[test]
1194        fn prop_adjoint_homomorphism(
1195            r1 in arb_so3(),
1196            r2 in arb_so3(),
1197            v in arb_so3_algebra()
1198        ) {
1199            // Compute Ad_{R₁∘R₂}(v)
1200            let r_composed = r1.compose(&r2);
1201            let left = r_composed.adjoint_action(&v);
1202
1203            // Compute Ad_{R₁}(Ad_{R₂}(v))
1204            let ad_r2_v = r2.adjoint_action(&v);
1205            let right = r1.adjoint_action(&ad_r2_v);
1206
1207            // They should be equal
1208            let diff = left.add(&right.scale(-1.0));
1209            prop_assert!(
1210                diff.norm() < 1e-7,
1211                "Adjoint homomorphism failed: Ad_{{R₁∘R₂}}(v) != Ad_{{R₁}}(Ad_{{R₂}}(v)), diff norm = {}",
1212                diff.norm()
1213            );
1214        }
1215
1216        /// **Adjoint Representation: Identity Action**
1217        ///
1218        /// The identity element acts trivially on the Lie algebra:
1219        /// - Ad_I(v) = v for all v ∈ so(3)
1220        #[test]
1221        fn prop_adjoint_identity(v in arb_so3_algebra()) {
1222            let e = SO3::identity();
1223            let result = e.adjoint_action(&v);
1224
1225            let diff = result.add(&v.scale(-1.0));
1226            prop_assert!(
1227                diff.norm() < 1e-10,
1228                "Identity action failed: Ad_I(v) != v, diff norm = {}",
1229                diff.norm()
1230            );
1231        }
1232
1233        /// **Adjoint Representation: Lie Bracket Preservation**
1234        ///
1235        /// The adjoint representation preserves the Lie bracket:
1236        /// - Ad_R([v,w]) = [Ad_R(v), Ad_R(w)]
1237        ///
1238        /// This is a critical property that ensures the adjoint action
1239        /// is a Lie algebra automorphism.
1240        #[test]
1241        fn prop_adjoint_bracket_preservation(
1242            r in arb_so3(),
1243            v in arb_so3_algebra(),
1244            w in arb_so3_algebra()
1245        ) {
1246            use crate::traits::LieAlgebra;
1247
1248            // Compute Ad_R([v,w])
1249            let bracket_vw = v.bracket(&w);
1250            let left = r.adjoint_action(&bracket_vw);
1251
1252            // Compute [Ad_R(v), Ad_R(w)]
1253            let ad_v = r.adjoint_action(&v);
1254            let ad_w = r.adjoint_action(&w);
1255            let right = ad_v.bracket(&ad_w);
1256
1257            // They should be equal
1258            let diff = left.add(&right.scale(-1.0));
1259            prop_assert!(
1260                diff.norm() < 1e-6,
1261                "Bracket preservation failed: Ad_R([v,w]) != [Ad_R(v), Ad_R(w)], diff norm = {}",
1262                diff.norm()
1263            );
1264        }
1265
1266        /// **Adjoint Representation: Inverse Property**
1267        ///
1268        /// The inverse of an element acts as the inverse transformation:
1269        /// - Ad_{R⁻¹}(Ad_R(v)) = v
1270        #[test]
1271        fn prop_adjoint_inverse(r in arb_so3(), v in arb_so3_algebra()) {
1272            // Apply Ad_R then Ad_{R⁻¹}
1273            let ad_r_v = r.adjoint_action(&v);
1274            let r_inv = r.inverse();
1275            let result = r_inv.adjoint_action(&ad_r_v);
1276
1277            // Should recover v
1278            let diff = result.add(&v.scale(-1.0));
1279            prop_assert!(
1280                diff.norm() < 1e-7,
1281                "Inverse property failed: Ad_{{R⁻¹}}(Ad_R(v)) != v, diff norm = {}",
1282                diff.norm()
1283            );
1284        }
1285    }
1286}