1use crate::traits::{AntiHermitianByConstruction, LieAlgebra, LieGroup, TracelessByConstruction};
40use nalgebra::Matrix3;
41use std::fmt;
42use std::ops::{Add, Mul, MulAssign, Neg, Sub};
43
44#[derive(Clone, Copy, Debug, PartialEq)]
67pub struct So3Algebra(pub(crate) [f64; 3]);
68
69impl So3Algebra {
70 #[must_use]
75 pub fn new(components: [f64; 3]) -> Self {
76 Self(components)
77 }
78
79 #[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 #[inline]
221 fn bracket(&self, other: &Self) -> Self {
222 let v = self.0;
224 let w = other.0;
225
226 Self([
227 v[1] * w[2] - v[2] * w[1], v[2] * w[0] - v[0] * w[2], v[0] * w[1] - v[1] * w[0], ])
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#[derive(Clone, Debug, PartialEq)]
265pub struct SO3 {
266 pub(crate) matrix: Matrix3<f64>,
268}
269
270impl SO3 {
271 #[must_use]
273 pub fn matrix(&self) -> &Matrix3<f64> {
274 &self.matrix
275 }
276
277 #[must_use]
279 pub fn identity() -> Self {
280 Self {
281 matrix: Matrix3::identity(),
282 }
283 }
284
285 #[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 #[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 #[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 #[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 let n = [axis[0] / norm, axis[1] / norm, axis[2] / norm];
344
345 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 #[must_use]
366 pub fn trace(&self) -> f64 {
367 self.matrix.trace()
368 }
369
370 #[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 #[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 #[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 #[must_use]
407 pub fn inverse(&self) -> Self {
408 Self {
409 matrix: self.matrix.transpose(),
410 }
411 }
412
413 #[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 let cos_theta = cos_theta.clamp(-1.0, 1.0);
426
427 cos_theta.acos()
428 }
429
430 #[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 let interpolated = self.matrix * (1.0 - t) + other.matrix * t;
457
458 Self::gram_schmidt_orthogonalize(interpolated)
460 }
461
462 #[must_use]
473 pub fn gram_schmidt_orthogonalize(matrix: Matrix3<f64>) -> Self {
474 use nalgebra::Vector3;
475
476 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 let e0 = c0.normalize();
482 let e1 = (c1 - e0 * e0.dot(&c1)).normalize();
483 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 #[must_use]
498 pub fn renormalize(&self) -> Self {
499 Self::gram_schmidt_orthogonalize(self.matrix)
500 }
501
502 #[must_use]
508 pub fn geodesic_distance(&self, other: &Self) -> f64 {
509 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
567impl 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 self.inverse()
612 }
613
614 fn adjoint_action(&self, algebra_element: &So3Algebra) -> So3Algebra {
615 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 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 let trace = self.matrix.trace();
658
659 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 return Ok(So3Algebra::zero());
669 }
670
671 if (theta - std::f64::consts::PI).abs() < SINGULARITY_THRESHOLD {
672 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 let sin_theta = theta.sin();
695
696 if sin_theta.abs() < 1e-10 {
697 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 Ok(So3Algebra([theta * nx, theta * ny, theta * nz]))
710 }
711}
712
713use crate::traits::{Compact, SemiSimple, Simple};
718
719impl Compact for SO3 {}
724
725impl Simple for SO3 {}
729
730impl SemiSimple for SO3 {}
732
733impl TracelessByConstruction for So3Algebra {}
742
743impl 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 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 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 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 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 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 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 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 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 let det = fixed.matrix.determinant();
896 assert_relative_eq!(det, 1.0, epsilon = 1e-10);
897 }
898
899 #[test]
900 fn test_renormalize() {
901 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 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 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 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 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 #[cfg(feature = "proptest")]
1021 use proptest::prelude::*;
1022
1023 #[cfg(feature = "proptest")]
1030 fn arb_so3() -> impl Strategy<Value = SO3> {
1031 use std::f64::consts::TAU;
1032
1033 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 #[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 #[test]
1069 fn prop_identity_axiom(r in arb_so3()) {
1070 let e = SO3::identity();
1071
1072 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 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 #[test]
1100 fn prop_inverse_axiom(r in arb_so3()) {
1101 let r_inv = r.inverse();
1102
1103 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 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 #[test]
1130 fn prop_associativity(r1 in arb_so3(), r2 in arb_so3(), r3 in arb_so3()) {
1131 let left_assoc = r1.compose(&r2).compose(&r3);
1133
1134 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 #[test]
1149 fn prop_inverse_continuity(r in arb_so3()) {
1150 let epsilon = 0.01;
1152 let perturbation = SO3::rotation_x(epsilon);
1153 let r_perturbed = r.compose(&perturbation);
1154
1155 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 #[test]
1170 fn prop_orthogonality_preserved(r1 in arb_so3(), r2 in arb_so3()) {
1171 let product = r1.compose(&r2);
1173 prop_assert!(
1174 product.verify_orthogonality(1e-10),
1175 "Composition violated orthogonality"
1176 );
1177
1178 let inv = r1.inverse();
1180 prop_assert!(
1181 inv.verify_orthogonality(1e-10),
1182 "Inverse violated orthogonality"
1183 );
1184 }
1185
1186 #[test]
1194 fn prop_adjoint_homomorphism(
1195 r1 in arb_so3(),
1196 r2 in arb_so3(),
1197 v in arb_so3_algebra()
1198 ) {
1199 let r_composed = r1.compose(&r2);
1201 let left = r_composed.adjoint_action(&v);
1202
1203 let ad_r2_v = r2.adjoint_action(&v);
1205 let right = r1.adjoint_action(&ad_r2_v);
1206
1207 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 #[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 #[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 let bracket_vw = v.bracket(&w);
1250 let left = r.adjoint_action(&bracket_vw);
1251
1252 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 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 #[test]
1271 fn prop_adjoint_inverse(r in arb_so3(), v in arb_so3_algebra()) {
1272 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 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}