IMP  2.1.0
The Integrative Modeling Platform
Rotation3D.h
Go to the documentation of this file.
1 /**
2  * \file IMP/algebra/Rotation3D.h \brief Simple 3D rotation class.
3  *
4  * Copyright 2007-2013 IMP Inventors. All rights reserved.
5  *
6  */
7 
8 #ifndef IMPALGEBRA_ROTATION_3D_H
9 #define IMPALGEBRA_ROTATION_3D_H
10 
11 #include <IMP/algebra/algebra_config.h>
12 #include "Vector3D.h"
13 #include "utility.h"
14 #include "constants.h"
15 #include "GeometricPrimitiveD.h"
16 
17 #include <IMP/base/log.h>
18 #include <cmath>
19 #include <iostream>
20 #include <algorithm>
21 
22 IMPALGEBRA_BEGIN_NAMESPACE
23 
24 #if !defined(IMP_DOXYGEN) && !defined(SWIG)
25 class Rotation3D;
26 Rotation3D compose(const Rotation3D &a, const Rotation3D &b);
27 #endif
28 
29 //! 3D rotation class.
30 /** Rotations are currently represented using quaternions and a cached
31  copy of the rotation matrix. The quaternion allows for fast and
32  stable composition and the cached rotation matrix means that
33  rotations are performed quickly. See
34  http://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation for
35  a comparison of different implementations of rotations.
36 
37  Currently the rotation can be initialized from either:
38  - XYZ Euler angles
39  - Rotation Matrix
40  - Quaternion
41  - angle/axis representation
42 
43  \geometry
44 */
45 class IMPALGEBRAEXPORT Rotation3D : public GeometricPrimitiveD<3> {
46  VectorD<4> v_;
47  mutable bool has_cache_;
48  mutable Vector3D matrix_[3];
49  IMP_NO_SWIG(friend Rotation3D compose(const Rotation3D &a,
50  const Rotation3D &b));
51  void fill_cache() const {
52  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
53  "Attempting to apply uninitialized rotation");
54  has_cache_ = true;
55  double v0s = get_squared(v_[0]);
56  double v1s = get_squared(v_[1]);
57  double v2s = get_squared(v_[2]);
58  double v3s = get_squared(v_[3]);
59  double v12 = v_[1] * v_[2];
60  double v01 = v_[0] * v_[1];
61  double v02 = v_[0] * v_[2];
62  double v23 = v_[2] * v_[3];
63  double v03 = v_[0] * v_[3];
64  double v13 = v_[1] * v_[3];
65  matrix_[0] =
66  Vector3D(v0s + v1s - v2s - v3s, 2 * (v12 - v03), 2 * (v13 + v02));
67  matrix_[1] =
68  Vector3D(2 * (v12 + v03), v0s - v1s + v2s - v3s, 2 * (v23 - v01));
69  matrix_[2] =
70  Vector3D(2 * (v13 - v02), 2 * (v23 + v01), v0s - v1s - v2s + v3s);
71  }
72 
73  public:
74  IMP_CXX11_DEFAULT_COPY_CONSTRUCTOR(Rotation3D);
75  //! Create a rotation from an unnormalized vector 4
77  : v_(v.get_unit_vector()), has_cache_(false) {}
78 
79  //! Create an invalid rotation
80  Rotation3D() : v_(0, 0, 0, 0), has_cache_(false) {}
81  //! Create a rotation from a quaternion
82  /** \throw base::ValueException if the rotation is not a unit vector.
83  */
84  Rotation3D(double a, double b, double c, double d)
85  : v_(a, b, c, d), has_cache_(false) {
87  v_.get_squared_magnitude(), 1.0,
88  "Attempting to construct a rotation from a "
89  << " non-quaternion value. The coefficient vector"
90  << " must have a length of 1. Got: " << a << " " << b << " " << c
91  << " " << d << " gives " << v_.get_squared_magnitude());
92  if (a < 0) {
93  // make them canonical
94  v_ = -v_;
95  }
96  }
97  ~Rotation3D();
98 
99 #ifndef IMP_DOXYGEN
100  Vector3D get_rotated_no_cache(const Vector3D &o) const {
101  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
102  "Attempting to access uninitialized rotation");
103  return Vector3D(
104  (v_[0] * v_[0] + v_[1] * v_[1] - v_[2] * v_[2] - v_[3] * v_[3]) * o[0] +
105  2 * (v_[1] * v_[2] - v_[0] * v_[3]) * o[1] +
106  2 * (v_[1] * v_[3] + v_[0] * v_[2]) * o[2],
107  2 * (v_[1] * v_[2] + v_[0] * v_[3]) * o[0] +
108  (v_[0] * v_[0] - v_[1] * v_[1] + v_[2] * v_[2] - v_[3] * v_[3]) *
109  o[1] +
110  2 * (v_[2] * v_[3] - v_[0] * v_[1]) * o[2],
111  2 * (v_[1] * v_[3] - v_[0] * v_[2]) * o[0] +
112  2 * (v_[2] * v_[3] + v_[0] * v_[1]) * o[1] +
113  (v_[0] * v_[0] - v_[1] * v_[1] - v_[2] * v_[2] + v_[3] * v_[3]) *
114  o[2]);
115  }
116 
117  //! Gets only the requested rotation coordinate of the vector
118  double get_rotated_one_coordinate_no_cache(const Vector3D &o,
119  unsigned int coord) const {
120  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
121  "Attempting to apply uninitialized rotation");
122  switch (coord) {
123  case 0:
124  return (v_[0] * v_[0] + v_[1] * v_[1] - v_[2] * v_[2] - v_[3] * v_[3]) *
125  o[0] +
126  2 * (v_[1] * v_[2] - v_[0] * v_[3]) * o[1] +
127  2 * (v_[1] * v_[3] + v_[0] * v_[2]) * o[2];
128  break;
129  case 1:
130  return 2 * (v_[1] * v_[2] + v_[0] * v_[3]) * o[0] +
131  (v_[0] * v_[0] - v_[1] * v_[1] + v_[2] * v_[2] - v_[3] * v_[3]) *
132  o[1] +
133  2 * (v_[2] * v_[3] - v_[0] * v_[1]) * o[2];
134 
135  break;
136  case 2:
137  return 2 * (v_[1] * v_[3] - v_[0] * v_[2]) * o[0] +
138  2 * (v_[2] * v_[3] + v_[0] * v_[1]) * o[1] +
139  (v_[0] * v_[0] - v_[1] * v_[1] - v_[2] * v_[2] + v_[3] * v_[3]) *
140  o[2];
141  break;
142  default:
143  IMP_THROW("Out of range coordinate " << coord, base::IndexException);
144  }
145  }
146 #endif
147  //! Rotate a vector around the origin
148  Vector3D get_rotated(const Vector3D &o) const {
149  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
150  "Attempting to apply uninitialized rotation");
151  if (!has_cache_) fill_cache();
152  return Vector3D(o * matrix_[0], o * matrix_[1], o * matrix_[2]);
153  }
154 
155  //! Gets only the requested rotation coordinate of the vector
157  unsigned int coord) const {
158  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
159  "Attempting to apply uninitialized rotation");
160  if (!has_cache_) fill_cache();
161  return o * matrix_[coord];
162  }
163 
164  //! Rotate a vector around the origin
165  Vector3D operator*(const Vector3D &v) const { return get_rotated(v); }
166  Vector3D get_rotation_matrix_row(int i) const {
167  IMP_USAGE_CHECK((i >= 0) && (i <= 2), "row index out of range");
168  if (!has_cache_) fill_cache();
169  return matrix_[i];
170  }
171  IMP_SHOWABLE_INLINE(Rotation3D, {
172  out << v_[0] << " " << v_[1] << " " << v_[2] << " " << v_[3];
173  });
174 
175  //! Return the rotation which undoes this rotation.
176  inline Rotation3D get_inverse() const {
177  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
178  "Attempting to invert uninitialized rotation");
179  Rotation3D ret(v_[0], -v_[1], -v_[2], -v_[3]);
180  return ret;
181  }
182 
183  //! Return the quaternion so that it can be stored
184  /** Note that there is no guarantee on which of the two
185  equivalent quaternions is returned.
186  */
187  const Vector4D &get_quaternion() const {
188  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
189  "Attempting to access uninitialized rotation");
190  return v_;
191  }
192  //! multiply two rotations
193  Rotation3D operator*(const Rotation3D &q) const {
194  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
195  "Attempting to compose uninitialized rotation");
196  return compose(*this, q);
197  }
198 
199  //! Compute the rotation which when composed with r gives this
200  Rotation3D operator/(const Rotation3D &r) const {
201  IMP_USAGE_CHECK(v_.get_squared_magnitude() > 0,
202  "Attempting to compose uninitialized rotation");
203  return compose(*this, r.get_inverse());
204  }
205 
206  const Rotation3D &operator/=(const Rotation3D &r) {
207  *this = *this / r;
208  return *this;
209  }
210 
211  /** \brief Return the derivative of the position x with respect to
212  internal variable i. */
213  const Vector3D get_derivative(const Vector3D &o, unsigned int i) const;
214 };
215 
217 
218 //! Return a rotation that does not do anything
219 /** See Rotation3D */
220 inline Rotation3D get_identity_rotation_3d() { return Rotation3D(1, 0, 0, 0); }
221 
222 //! Return a distance between the two rotations
223 /** The distance runs between 0 and 1. More precisely,
224  the distance returned is distance between the two
225  quaternion vectors properly normalized, divided
226  by sqrt(2).
227 
228  A vector with distance d from the unit vector
229  represents a rotation of
230 
231  \f$ \theta= \cos^{-1}\left(1-\sqrt{2}d\right)\f$
232  See Rotation3D
233 */
234 inline double get_distance(const Rotation3D &r0, const Rotation3D &r1) {
235  double dot =
236  (r0.get_quaternion() - r1.get_quaternion()).get_squared_magnitude();
237  double odot =
238  (r0.get_quaternion() + r1.get_quaternion()).get_squared_magnitude();
239  double ans = std::min(dot, odot);
240  const double s2 = std::sqrt(2.0);
241  double ret = ans / s2;
242  return std::max(std::min(ret, 1.0), 0.0);
243 }
244 
245 //! Generate a Rotation3D object from a rotation around an axis
246 //! that is assumed to be normalized
247 /**
248  \param[in] axis_norm the normalized rotation axis passing through (0,0,0)
249  \param[in] angle the rotation angle in radians in the
250  clockwise direction
251  \note http://en.wikipedia.org/wiki/Rotation_matrix
252  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
253  angleToQuaternion/index.htm
254  See Rotation3D
255 */
256 IMPALGEBRAEXPORT Rotation3D
257  get_rotation_about_normalized_axis(const Vector3D &axis_norm, double angle);
258 
259 //! Generate a Rotation3D object from a rotation around an axis
260 /**
261  \param[in] axis the rotation axis passes through (0,0,0)
262  \param[in] angle the rotation angle in radians in the
263  clockwise direction
264  \note http://en.wikipedia.org/wiki/Rotation_matrix
265  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
266  angleToQuaternion/index.htm
267  See Rotation3D
268 */
269 IMPALGEBRAEXPORT Rotation3D
270  get_rotation_about_axis(const Vector3D &axis, double angle);
271 
272 //! Create a rotation from the first vector to the second one.
273 /** See Rotation3D
274  */
275 IMPALGEBRAEXPORT Rotation3D
277 
278 //! Generate a Rotation3D object from a rotation matrix
279 /**
280  See Rotation3D
281 */
282 IMPALGEBRAEXPORT Rotation3D
283  get_rotation_from_matrix(double m00, double m01, double m02, double m10,
284  double m11, double m12, double m20, double m21,
285  double m22);
286 
287 //! Pick a rotation at random from all possible rotations
288 /** See Rotation3D */
289 IMPALGEBRAEXPORT Rotation3D get_random_rotation_3d();
290 
291 //! Pick a rotation at random near the provided one
292 /** This method generates a rotation that is within the provided
293  distance of center.
294  \param[in] center The center of the rotational volume
295  \param[in] distance See
296  get_distance(const Rotation3D&,const Rotation3D&)
297  for a full definition.
298 
299  \note The cost of this operation increases as distance goes to 0.
300 
301  See Rotation3D
302 */
303 IMPALGEBRAEXPORT Rotation3D
304  get_random_rotation_3d(const Rotation3D &center, double distance);
305 
306 //! Cover the space of rotations evenly
307 /** If you care about the distance between samples instead of the number
308  of samples, the "surface area" of the set of rotations is pi^2. If
309  you allocate each sample a volume of 4/3 pi d^3 (to space them d apart),
310  Then you want 3/4 pi/d^3 points.
311 
312  Creates at least num_points rotations.
313 */
314 IMPALGEBRAEXPORT Rotation3Ds
315  get_uniform_cover_rotations_3d(unsigned int num_points);
316 
317 //! Generates a nondegenerate set of Euler angles with a delta resolution
318 /**
319 \param[in] delta sample every delta angles in radians.
320  */
322  double delta);
323 
324 //! Compute a rotatation from an unnormalized quaternion
325 /** See Rotation3D */
327  VectorD<4> uv = v.get_unit_vector();
328  return Rotation3D(uv[0], uv[1], uv[2], uv[3]);
329 }
330 
331 /** See Rotation3D
332  */
333 inline Rotation3D compose(const Rotation3D &a, const Rotation3D &b) {
334  return Rotation3D(a.v_[0] * b.v_[0] - a.v_[1] * b.v_[1] - a.v_[2] * b.v_[2] -
335  a.v_[3] * b.v_[3],
336  a.v_[0] * b.v_[1] + a.v_[1] * b.v_[0] + a.v_[2] * b.v_[3] -
337  a.v_[3] * b.v_[2],
338  a.v_[0] * b.v_[2] - a.v_[1] * b.v_[3] + a.v_[2] * b.v_[0] +
339  a.v_[3] * b.v_[1],
340  a.v_[0] * b.v_[3] + a.v_[1] * b.v_[2] - a.v_[2] * b.v_[1] +
341  a.v_[3] * b.v_[0]);
342 }
343 
344 /** \name Euler Angles
345  There are many conventions for how to define Euler angles, based on choices
346  of which of the x,y,z axis to use in what order and whether the rotation
347  axis is in the body frame (and hence affected by previous rotations) or in
348  in a fixed frame. See
349  http://en.wikipedia.org/wiki/Conversion_between_quaternions_and_Euler_angles
350  for a general description.
351 
352  - All Euler angles are specified in radians.
353  - Only x-y-z order is currently supported.
354  - To convert Euler angles in a different order, one can compose a Rotation3D
355  from three rotations using get_rotation_about_axis function.
356  @{
357 */
358 
359 //! Initialize a rotation in x-y-z order from three angles
360 /** \param[in] xr Rotation around the X axis in radians
361  \param[in] yr Rotation around the Y axis in radians
362  \param[in] zr Rotation around the Z axis in radians
363  \note The three rotations are represented in the original (fixed)
364  coordinate frame.
365  See Rotation3D
366  See FixedXYZ
367 */
368 IMPALGEBRAEXPORT Rotation3D
369  get_rotation_from_fixed_xyz(double xr, double yr, double zr);
370 
371 //! Initialize a rotation from euler angles
372 /**
373  \param[in] phi Rotation around the Z axis in radians
374  \param[in] theta Rotation around the X axis in radians
375  \param[in] psi Rotation around the Z axis in radians
376  \note The first rotation is by an angle phi about the z-axis.
377  The second rotation is by an angle theta in [0,pi] about the
378  former x-axis , and the third rotation is by an angle psi
379  about the former z-axis.
380  See Rotation3D
381 */
382 IMPALGEBRAEXPORT Rotation3D
383  get_rotation_from_fixed_zxz(double phi, double theta, double psi);
384 
385 //! Generate a rotation object from Euler Angles
386 /** \note The first rotation is by an angle about the z-axis.
387  The second rotation is by an angle about the new y-axis.
388  The third rotation is by an angle about the new z-axis.
389  \param[in] Rot First Euler angle (radians) defining the rotation (Z axis)
390  \param[in] Tilt Second Euler angle (radians) defining the rotation (Y axis)
391  \param[in] Psi Third Euler angle (radians) defining the rotation (Z axis)
392  See Rotation3D
393 */
394 IMPALGEBRAEXPORT Rotation3D
395  get_rotation_from_fixed_zyz(double Rot, double Tilt, double Psi);
396 
397 //! A simple class for returning XYZ Euler angles
398 class FixedXYZ : public GeometricPrimitiveD<3> {
399  double v_[3];
400 
401  public:
402  FixedXYZ() {}
403  FixedXYZ(double x, double y, double z) {
404  v_[0] = x;
405  v_[1] = y;
406  v_[2] = z;
407  }
408  double get_x() const { return v_[0]; }
409  double get_y() const { return v_[1]; }
410  double get_z() const { return v_[2]; }
412  { out << v_[0] << " " << v_[1] << " " << v_[2]; });
413 };
414 
416 
417 //! The inverse of rotation_from_fixed_xyz()
418 /**
419  \see rotation_from_fixed_xyz()
420  See Rotation3D
421  See FixesXYZ
422 */
423 IMPALGEBRAEXPORT FixedXYZ get_fixed_xyz_from_rotation(const Rotation3D &r);
424 
425 /** @}*/
426 
427 //! Interpolate between two rotations
428 /** It f ==0, return b, if f==1 return a.
429  See Rotation3D*/
431  double f) {
432  VectorD<4> bq = b.get_quaternion(), aq = a.get_quaternion();
433  if (bq * aq < 0) bq = -bq;
434  return f * aq + (1 - f) * bq;
435 }
436 
437 /** Return the rotation which takes the native x and y axes to the
438  given x and y axes.
439  The two axis must be perpendicular unit vectors.
440 */
441 IMPALGEBRAEXPORT Rotation3D
442  get_rotation_from_x_y_axes(const Vector3D &x, const Vector3D &y);
443 
444 //! Decompose a Rotation3D object into a rotation around an axis
445 /**
446  \note http://en.wikipedia.org/wiki/Rotation_matrix
447  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
448  angleToQuaternion/index.htm
449  See Rotation3D
450 
451  @return axis direction and rotation about the axis in Radians
452 */
453 IMPALGEBRAEXPORT std::pair<Vector3D, double> get_axis_and_angle(
454  const Rotation3D &rot);
455 
456 typedef std::pair<Vector3D, double> AxisAnglePair;
457 #ifndef IMP_DOXYGEN
458 typedef base::Vector<AxisAnglePair> AxisAnglePairs;
459 #endif
460 
461 IMPALGEBRA_END_NAMESPACE
462 #endif /* IMPALGEBRA_ROTATION_3D_H */
Vector3D get_rotated(const Vector3D &o) const
Rotate a vector around the origin.
Definition: Rotation3D.h:148
An exception for a request for an invalid member of a container.
Rotation3D get_rotation_about_normalized_axis(const Vector3D &axis_norm, double angle)
Basic types used by IMP.
Rotation3D get_rotation_from_x_y_axes(const Vector3D &x, const Vector3D &y)
Rotation3D get_rotation_from_fixed_xyz(double xr, double yr, double zr)
Initialize a rotation in x-y-z order from three angles.
Rotation3D()
Create an invalid rotation.
Definition: Rotation3D.h:80
Rotation2D compose(const Rotation2D &a, const Rotation2D &b)
compose two rotations a and b
Definition: Rotation2D.h:105
algebra::Rotation3Ds get_uniformly_sampled_rotations(double delta)
Generates a nondegenerate set of Euler angles with a delta resolution.
Rotation3D operator/(const Rotation3D &r) const
Compute the rotation which when composed with r gives this.
Definition: Rotation3D.h:200
A simple class for returning XYZ Euler angles.
Definition: Rotation3D.h:398
#define IMP_VALUES(Name, PluralName)
Define the type for storing sets of values.
Vector3D operator*(const Vector3D &v) const
Rotate a vector around the origin.
Definition: Rotation3D.h:165
#define IMP_SHOWABLE_INLINE(Name, how_to_show)
Declare the methods needed by an object that can be printed.
double get_rotated_one_coordinate(const Vector3D &o, unsigned int coord) const
Gets only the requested rotation coordinate of the vector.
Definition: Rotation3D.h:156
Rotation3D(const VectorD< 4 > &v)
Create a rotation from an unnormalized vector 4.
Definition: Rotation3D.h:76
#define IMP_USAGE_CHECK_FLOAT_EQUAL(expra, exprb, message)
const Vector4D & get_quaternion() const
Return the quaternion so that it can be stored.
Definition: Rotation3D.h:187
#define IMP_USAGE_CHECK(expr, message)
A runtime test for incorrect usage of a class or method.
Rotation3D get_rotation_about_axis(const Vector3D &axis, double angle)
Generate a Rotation3D object from a rotation around an axis.
IMP::base::Vector< Rotation3D > Rotation3Ds
Definition: Rotation3D.h:216
Rotation3D get_random_rotation_3d()
Pick a rotation at random from all possible rotations.
double get_distance(const Plane3D &pln, const Vector3D &p)
Return the distance between a plane and a point in 3D.
Definition: Plane3D.h:68
Functions to deal with very common math operations.
FixedXYZ get_fixed_xyz_from_rotation(const Rotation3D &r)
The inverse of rotation_from_fixed_xyz()
Rotation3D get_inverse() const
Return the rotation which undoes this rotation.
Definition: Rotation3D.h:176
Rotation3D get_rotation_from_fixed_zyz(double Rot, double Tilt, double Psi)
Generate a rotation object from Euler Angles.
#define IMP_NO_SWIG(x)
Hide the line when SWIG is compiled or parses it.
std::pair< Vector3D, double > get_axis_and_angle(const Rotation3D &rot)
Decompose a Rotation3D object into a rotation around an axis.
Rotation3D get_interpolated(const Rotation3D &a, const Rotation3D &b, double f)
Interpolate between two rotations.
Definition: Rotation3D.h:430
Rotation3D get_rotation_from_vector4d(const VectorD< 4 > &v)
Compute a rotatation from an unnormalized quaternion.
Definition: Rotation3D.h:326
Rotation3D get_rotation_taking_first_to_second(const Vector3D &v1, const Vector3D &v2)
Create a rotation from the first vector to the second one.
Rotation3Ds get_uniform_cover_rotations_3d(unsigned int num_points)
Cover the space of rotations evenly.
Various useful constants.
Rotation3D operator*(const Rotation3D &q) const
multiply two rotations
Definition: Rotation3D.h:193
3D rotation class.
Definition: Rotation3D.h:45
Rotation3D get_rotation_from_matrix(double m00, double m01, double m02, double m10, double m11, double m12, double m20, double m21, double m22)
Generate a Rotation3D object from a rotation matrix.
Rotation3D(double a, double b, double c, double d)
Create a rotation from a quaternion.
Definition: Rotation3D.h:84
#define IMP_THROW(message, exception_name)
Throw an exception with a message.
Logging and error reporting support.
VectorD< 3 > Vector3D
Definition: VectorD.h:587
Rotation3D get_identity_rotation_3d()
Return a rotation that does not do anything.
Definition: Rotation3D.h:220
Simple 3D vector class.
VectorD get_unit_vector() const
Definition: VectorD.h:204
Rotation3D get_rotation_from_fixed_zxz(double phi, double theta, double psi)
Initialize a rotation from euler angles.