28 #ifndef TOON_INCLUDE_SO3_H 29 #define TOON_INCLUDE_SO3_H 31 #include <TooN/TooN.h> 32 #include <TooN/helpers.h> 37 template <
typename Precision>
class SO3;
38 template <
typename Precision>
class SE3;
39 template <
typename Precision>
class SIM3;
41 template<
class Precision>
inline std::istream & operator>>(std::istream &,
SO3<Precision> & );
42 template<
class Precision>
inline std::istream & operator>>(std::istream &,
SE3<Precision> & );
43 template<
class Precision>
inline std::istream & operator>>(std::istream &,
SIM3<Precision> & );
52 template <
typename Precision = DefaultPrecision>
55 friend std::istream&
operator>> <Precision> (std::istream& is,
SO3<Precision> & rhs);
56 friend std::istream&
operator>> <Precision> (std::istream& is,
SE3<Precision> & rhs);
57 friend std::istream&
operator>> <Precision> (std::istream& is,
SIM3<Precision> & rhs);
58 friend class SE3<Precision>;
59 friend class SIM3<Precision>;
62 SO3() : my_matrix(Identity) {}
65 template <
int S,
typename P,
typename A>
69 template <
int R,
int C,
typename P,
typename A>
77 template <
int S1,
int S2,
typename P1,
typename P2,
typename A1,
typename A2>
93 R1.T()[2] = R1.T()[0] ^ n;
94 my_matrix.T()[0] =
unit(b);
96 my_matrix.T()[2] = my_matrix.T()[0] ^ n;
97 my_matrix = my_matrix * R1.T();
102 template <
int R,
int C,
typename P,
typename A>
111 my_matrix[0] =
unit(my_matrix[0]);
112 my_matrix[1] -= my_matrix[0] * (my_matrix[0]*my_matrix[1]);
113 my_matrix[1] =
unit(my_matrix[1]);
114 my_matrix[2] -= my_matrix[0] * (my_matrix[0]*my_matrix[2]);
115 my_matrix[2] -= my_matrix[1] * (my_matrix[1]*my_matrix[2]);
116 my_matrix[2] =
unit(my_matrix[2]);
118 assert( (my_matrix[0] ^ my_matrix[1]) * my_matrix[2] > 0 );
133 template <
typename P>
154 result[(i+1)%3][(i+2)%3] = -1;
155 result[(i+2)%3][(i+1)%3] = 1;
160 template<
typename Base>
165 result[(i+1)%3] = - pos[(i+2)%3];
166 result[(i+2)%3] = pos[(i+1)%3];
174 template <
int S,
typename A>
181 template <
typename PA,
typename PB>
186 inline SO3(
const SO3& so3,
const Invert&) : my_matrix(so3.my_matrix.T()) {}
193 template <
typename Precision>
194 inline std::ostream& operator<< (std::ostream& os, const SO3<Precision>& rhs){
195 return os << rhs.get_matrix();
200 template <
typename Precision>
218 template <
typename Precision,
int S,
typename VP,
typename VA,
typename MA>
222 const Precision wx2 = (Precision)w[0]*w[0];
223 const Precision wy2 = (Precision)w[1]*w[1];
224 const Precision wz2 = (Precision)w[2]*w[2];
226 R[0][0] = 1.0 - B*(wy2 + wz2);
227 R[1][1] = 1.0 - B*(wx2 + wz2);
228 R[2][2] = 1.0 - B*(wx2 + wy2);
231 const Precision a = A*w[2];
232 const Precision b = B*(w[0]*w[1]);
237 const Precision a = A*w[1];
238 const Precision b = B*(w[0]*w[2]);
243 const Precision a = A*w[0];
244 const Precision b = B*(w[1]*w[2]);
252 template <
typename Precision>
253 template<
int S,
typename VP,
typename VA>
260 static const Precision one_6th = 1.0/6.0;
261 static const Precision one_20th = 1.0/20.0;
265 const Precision theta_sq = w*w;
266 const Precision theta =
sqrt(theta_sq);
270 if (theta_sq < 1e-8) {
271 A = 1.0 - one_6th * theta_sq;
274 if (theta_sq < 1e-6) {
275 B = 0.5 - 0.25 * one_6th * theta_sq;
276 A = 1.0 - theta_sq * one_6th*(1.0 - one_20th * theta_sq);
278 const Precision inv_theta = 1.0/theta;
279 A = sin(theta) * inv_theta;
280 B = (1 - cos(theta)) * (inv_theta * inv_theta);
287 template <
typename Precision>
292 const Precision cos_angle = (my_matrix[0][0] + my_matrix[1][1] + my_matrix[2][2] - 1.0) * 0.5;
293 result[0] = (my_matrix[2][1]-my_matrix[1][2])/2;
294 result[1] = (my_matrix[0][2]-my_matrix[2][0])/2;
295 result[2] = (my_matrix[1][0]-my_matrix[0][1])/2;
297 Precision sin_angle_abs =
sqrt(result*result);
298 if (cos_angle > M_SQRT1_2) {
299 if(sin_angle_abs > 0){
300 result *= asin(sin_angle_abs) / sin_angle_abs;
302 }
else if( cos_angle > -M_SQRT1_2) {
303 const Precision angle = acos(cos_angle);
304 result *= angle / sin_angle_abs;
307 const Precision angle = M_PI - asin(sin_angle_abs);
308 const Precision d0 = my_matrix[0][0] - cos_angle,
309 d1 = my_matrix[1][1] - cos_angle,
310 d2 = my_matrix[2][2] - cos_angle;
312 if(d0*d0 > d1*d1 && d0*d0 > d2*d2){
314 r2[1] = (my_matrix[1][0]+my_matrix[0][1])/2;
315 r2[2] = (my_matrix[0][2]+my_matrix[2][0])/2;
316 }
else if(d1*d1 > d2*d2) {
317 r2[0] = (my_matrix[1][0]+my_matrix[0][1])/2;
319 r2[2] = (my_matrix[2][1]+my_matrix[1][2])/2;
321 r2[0] = (my_matrix[0][2]+my_matrix[2][0])/2;
322 r2[1] = (my_matrix[2][1]+my_matrix[1][2])/2;
329 result = TooN::operator*(angle,r2);
336 template<
int S,
typename P,
typename PV,
typename A>
inline 343 template<
int S,
typename P,
typename PV,
typename A>
inline 350 template<
int R,
int C,
typename P,
typename PM,
typename A>
inline 357 template<
int R,
int C,
typename P,
typename PM,
typename A>
inline 362 #if 0 // will be moved to another header file ? 364 template <
class A>
inline 365 Vector<3> transform(
const SO3& pose,
const FixedVector<3,A>& x) {
return pose*x; }
367 template <
class A1,
class A2>
inline 368 Vector<3> transform(
const SO3& pose,
const FixedVector<3,A1>& x, FixedMatrix<3,3,A2>& J_x)
374 template <
class A1,
class A2,
class A3>
inline 375 Vector<3> transform(
const SO3& pose,
const FixedVector<3,A1>& x, FixedMatrix<3,3,A2>& J_x, FixedMatrix<3,3,A3>& J_pose)
379 J_pose[0][0] = J_pose[1][1] = J_pose[2][2] = 0;
380 J_pose[1][0] = -(J_pose[0][1] = cx[2]);
381 J_pose[0][2] = -(J_pose[2][0] = cx[1]);
382 J_pose[2][1] = -(J_pose[1][2] = cx[0]);
386 template <
class A1,
class A2,
class A3>
inline 387 Vector<2> project_transformed_point(
const SO3& pose,
const FixedVector<3,A1>& in_frame, FixedMatrix<2,3,A2>& J_x, FixedMatrix<2,3,A3>& J_pose)
389 const double z_inv = 1.0/in_frame[2];
390 const double x_z_inv = in_frame[0]*z_inv;
391 const double y_z_inv = in_frame[1]*z_inv;
392 const double cross = x_z_inv * y_z_inv;
393 J_pose[0][0] = -cross;
394 J_pose[0][1] = 1 + x_z_inv*x_z_inv;
395 J_pose[0][2] = -y_z_inv;
396 J_pose[1][0] = -1 - y_z_inv*y_z_inv;
397 J_pose[1][1] = cross;
398 J_pose[1][2] = x_z_inv;
401 J_x[0][0] = z_inv*(R[0][0] - x_z_inv * R[2][0]);
402 J_x[0][1] = z_inv*(R[0][1] - x_z_inv * R[2][1]);
403 J_x[0][2] = z_inv*(R[0][2] - x_z_inv * R[2][2]);
404 J_x[1][0] = z_inv*(R[1][0] - y_z_inv * R[2][0]);
405 J_x[1][1] = z_inv*(R[1][1] - y_z_inv * R[2][1]);
406 J_x[1][2] = z_inv*(R[1][2] - y_z_inv * R[2][2]);
408 return makeVector(x_z_inv, y_z_inv);
412 template <
class A1>
inline 413 Vector<2> transform_and_project(
const SO3& pose,
const FixedVector<3,A1>& x)
415 return project(transform(pose,x));
418 template <
class A1,
class A2,
class A3>
inline 419 Vector<2> transform_and_project(
const SO3& pose,
const FixedVector<3,A1>& x, FixedMatrix<2,3,A2>& J_x, FixedMatrix<2,3,A3>& J_pose)
421 return project_transformed_point(pose, transform(pose,x), J_x, J_pose);
Vector< 3, Precision > ln() const
Take the logarithm of the matrix, generating the corresponding vector in the Lie Algebra.
Definition: so3.h:288
Vector< 3, typename Internal::MultiplyType< PV, P >::type > operator*(const Vector< S, PV, A > &lhs, const SO3< P > &rhs)
Left-multiply by a Vector.
Definition: so3.h:344
Represent a three-dimensional Euclidean transformation (a rotation and a translation).
Definition: se3.h:48
SO3 & operator*=(const SO3< P > &rhs)
Right-multiply by another rotation matrix.
Definition: so3.h:134
void rodrigues_so3_exp(const Vector< S, VP, VA > &w, const Precision A, const Precision B, Matrix< 3, 3, Precision, MA > &R)
Compute a rotation exponential using the Rodrigues Formula.
Definition: so3.h:219
SO3 & operator=(const Matrix< R, C, P, A > &rhs)
Assignment operator from a general matrix.
Definition: so3.h:103
Vector< 3, typename Internal::MultiplyType< P, PV >::type > operator*(const SO3< P > &lhs, const Vector< S, PV, A > &rhs)
Right-multiply by a Vector.
Definition: so3.h:337
Vector< 3, Precision > adjoint(const Vector< S, Precision, A > &vect) const
Transfer a vector in the Lie Algebra from one co-ordinate frame to another such that for a matrix ...
Definition: so3.h:175
SO3()
Default constructor. Initialises the matrix to the identity (no rotation)
Definition: so3.h:62
Vector<(Size==Dynamic?Dynamic:Size-1)+0, Precision > project(const Vector< Size, Precision, Base > &v)
For a vector v of length i, return .
Definition: helpers.h:182
static SO3 exp(const Vector< S, VP, A > &vect)
Exponentiate a vector in the Lie algebra to generate a new SO3.
Pretty generic SFINAE introspection generator.
Definition: vec_test.cc:21
std::istream & operator>>(std::istream &is, SO3< Precision > &rhs)
Read from SO3 to a stream.
Definition: so3.h:201
SO3 inverse() const
Returns the inverse of this matrix (=the transpose, so this is a fast operation)
Definition: so3.h:130
Precision norm_sq(const Vector< Size, Precision, Base > &v)
Compute the norm of v.
Definition: helpers.h:106
static Vector< 3, Precision > generator_field(int i, const Vector< 3, Precision, Base > &pos)
Returns the i-th generator times pos.
Definition: so3.h:161
static Matrix< 3, 3, Precision > generator(int i)
Returns the i-th generator.
Definition: so3.h:152
SO3< typename Internal::MultiplyType< Precision, P >::type > operator*(const SO3< P > &rhs) const
Right-multiply by another rotation matrix.
Definition: so3.h:141
Matrix< 3, C, typename Internal::MultiplyType< P, PM >::type > operator*(const SO3< P > &lhs, const Matrix< R, C, PM, A > &rhs)
Right-multiply by a matrix.
Definition: so3.h:351
Class to represent a three-dimensional rotation matrix.
Definition: so3.h:37
SO3(const Vector< S1, P1, A1 > &a, const Vector< S2, P2, A2 > &b)
creates an SO3 as a rotation that takes Vector a into the direction of Vector b with the rotation axi...
Definition: so3.h:78
Matrix< R, 3, typename Internal::MultiplyType< PM, P >::type > operator*(const Matrix< R, C, PM, A > &lhs, const SO3< P > &rhs)
Left-multiply by a matrix.
Definition: so3.h:358
Represent a three-dimensional similarity transformation (a rotation, a scale factor and a translation...
Definition: sim3.h:48
Matrix< R, C, P > sqrt(const Matrix< R, C, P, B > &m)
computes a matrix square root of a matrix m by the product form of the Denman and Beavers iteration a...
Definition: helpers.h:350
Definition: size_mismatch.hh:103
void coerce()
Modifies the matrix to make sure it is a valid rotation matrix.
Definition: so3.h:110
SO3(const Matrix< R, C, P, A > &rhs)
Construct from a rotation matrix.
Definition: so3.h:70
const Matrix< 3, 3, Precision > & get_matrix() const
Returns the SO3 as a Matrix<3>
Definition: so3.h:146
Vector< Size, Precision > unit(const Vector< Size, Precision, Base > &v)
Compute a the unit vector .
Definition: helpers.h:153
SO3(const Vector< S, P, A > &v)
Construct from the axis of rotation (and angle given by the magnitude).
Definition: so3.h:66