IMP logo
IMP Reference Guide  2.22.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-2022 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 #include <cereal/access.hpp>
17 #include <Eigen/Dense>
18 
19 #include <IMP/log.h>
20 #include <cmath>
21 #include <iostream>
22 #include <algorithm>
23 
24 IMPALGEBRA_BEGIN_NAMESPACE
25 
26 typedef Vector4D Rotation3DAdjoint;
27 typedef std::pair<Vector3D, Rotation3DAdjoint> RotatedVector3DAdjoint;
28 typedef std::pair<Rotation3DAdjoint, Rotation3DAdjoint>
29  ComposeRotation3DAdjoint;
30 
31 #if !defined(IMP_DOXYGEN) && !defined(SWIG)
32 class Rotation3D;
33 Rotation3D compose(const Rotation3D &a, const Rotation3D &b);
34 #endif
35 
36 //! 3D rotation class.
37 /** Rotations are currently represented using quaternions and a cached
38  copy of the rotation matrix. The quaternion allows for fast and
39  stable composition and the cached rotation matrix means that
40  rotations are performed quickly. See
41  https://en.wikipedia.org/wiki/Quaternions_and_spatial_rotation for
42  a comparison of different implementations of rotations.
43 
44  Currently the rotation can be initialized from either:
45  - XYZ Euler angles
46  - Rotation Matrix
47  - Quaternion
48  - angle/axis representation
49 
50  \geometry
51 */
52 class IMPALGEBRAEXPORT Rotation3D : public GeometricPrimitiveD<3> {
53  VectorD<4> v_;
54  mutable bool has_cache_;
55  mutable Vector3D matrix_[3];
56 
57  friend class cereal::access;
58 
59  template<class Archive>
60  void save(Archive &ar) const {
61  ar(v_);
62  }
63 
64  template<class Archive>
65  void load(Archive &ar) {
66  ar(v_);
67  has_cache_ = false;
68  }
69 
70  IMP_NO_SWIG(friend Rotation3D compose(const Rotation3D &a,
71  const Rotation3D &b));
72  void fill_cache() const {
73  IMP_USAGE_CHECK(get_is_valid(),
74  "Attempting to apply uninitialized rotation");
75  has_cache_ = true;
76  double v0s = get_squared(v_[0]);
77  double v1s = get_squared(v_[1]);
78  double v2s = get_squared(v_[2]);
79  double v3s = get_squared(v_[3]);
80  double v12 = v_[1] * v_[2];
81  double v01 = v_[0] * v_[1];
82  double v02 = v_[0] * v_[2];
83  double v23 = v_[2] * v_[3];
84  double v03 = v_[0] * v_[3];
85  double v13 = v_[1] * v_[3];
86  matrix_[0] =
87  Vector3D(v0s + v1s - v2s - v3s, 2 * (v12 - v03), 2 * (v13 + v02));
88  matrix_[1] =
89  Vector3D(2 * (v12 + v03), v0s - v1s + v2s - v3s, 2 * (v23 - v01));
90  matrix_[2] =
91  Vector3D(2 * (v13 - v02), 2 * (v23 + v01), v0s - v1s - v2s + v3s);
92  }
93 
94  // When copying preserve the cached rotation matrix, if present
95  void copy_cache(const Rotation3D &rot) {
96  has_cache_ = rot.has_cache_;
97  if (has_cache_) {
98  matrix_[0] = rot.matrix_[0];
99  matrix_[1] = rot.matrix_[1];
100  matrix_[2] = rot.matrix_[2];
101  }
102  }
103 
104  public:
105  //! Rotation3D copy constructor
106  Rotation3D(const Rotation3D &rot) : v_(rot.v_) {
107  copy_cache(rot);
108  }
109 
110  Rotation3D &operator=(const Rotation3D &rot) {
111  v_ = rot.v_;
112  copy_cache(rot);
113  return *this;
114  }
115 
116  //! Create a rotation from a vector of 4 quaternion coefficients.
117  //! @note: use assume_normalized with care - inputting an unnormalized
118  //! vector would result in unexpected results if it is true
119  explicit Rotation3D(const VectorD<4> &v,
120  bool assume_normalized=false)
121  : v_(assume_normalized ? v : v.get_unit_vector()),
122  has_cache_(false) {}
123 
124  //! Create an invalid rotation
125  Rotation3D() : v_(0, 0, 0, 0), has_cache_(false) {}
126  //! Create a rotation from a quaternion
127  /** \throw ValueException if the rotation is not a unit vector.
128  */
129  Rotation3D(double a, double b, double c, double d)
130  : v_(a, b, c, d), has_cache_(false) {
132  v_.get_squared_magnitude(), 1.0,
133  "Attempting to construct a rotation from a "
134  << " non-quaternion value. The coefficient vector"
135  << " must have a length of 1. Got: " << a << " " << b << " " << c
136  << " " << d << " gives " << v_.get_squared_magnitude());
137  if (a < 0) {
138  // make them canonical
139  v_ = -v_;
140  }
141  }
142  ~Rotation3D(){}
143 
144 #ifndef IMP_DOXYGEN
145  Vector3D get_rotated_no_cache(const Vector3D &o) const {
146  IMP_USAGE_CHECK(get_is_valid(),
147  "Attempting to access uninitialized rotation");
148  return Vector3D(
149  (v_[0] * v_[0] + v_[1] * v_[1] - v_[2] * v_[2] - v_[3] * v_[3]) * o[0] +
150  2 * (v_[1] * v_[2] - v_[0] * v_[3]) * o[1] +
151  2 * (v_[1] * v_[3] + v_[0] * v_[2]) * o[2],
152  2 * (v_[1] * v_[2] + v_[0] * v_[3]) * o[0] +
153  (v_[0] * v_[0] - v_[1] * v_[1] + v_[2] * v_[2] - v_[3] * v_[3]) *
154  o[1] +
155  2 * (v_[2] * v_[3] - v_[0] * v_[1]) * o[2],
156  2 * (v_[1] * v_[3] - v_[0] * v_[2]) * o[0] +
157  2 * (v_[2] * v_[3] + v_[0] * v_[1]) * o[1] +
158  (v_[0] * v_[0] - v_[1] * v_[1] - v_[2] * v_[2] + v_[3] * v_[3]) *
159  o[2]);
160  }
161 
162  //! Get only the requested rotation coordinate of the vector
163  double get_rotated_one_coordinate_no_cache(const Vector3D &o,
164  unsigned int coord) const {
165  IMP_USAGE_CHECK(get_is_valid(),
166  "Attempting to apply uninitialized rotation");
167  switch (coord) {
168  case 0:
169  return (v_[0] * v_[0] + v_[1] * v_[1] - v_[2] * v_[2] - v_[3] * v_[3]) *
170  o[0] +
171  2 * (v_[1] * v_[2] - v_[0] * v_[3]) * o[1] +
172  2 * (v_[1] * v_[3] + v_[0] * v_[2]) * o[2];
173  break;
174  case 1:
175  return 2 * (v_[1] * v_[2] + v_[0] * v_[3]) * o[0] +
176  (v_[0] * v_[0] - v_[1] * v_[1] + v_[2] * v_[2] - v_[3] * v_[3]) *
177  o[1] +
178  2 * (v_[2] * v_[3] - v_[0] * v_[1]) * o[2];
179 
180  break;
181  case 2:
182  return 2 * (v_[1] * v_[3] - v_[0] * v_[2]) * o[0] +
183  2 * (v_[2] * v_[3] + v_[0] * v_[1]) * o[1] +
184  (v_[0] * v_[0] - v_[1] * v_[1] - v_[2] * v_[2] + v_[3] * v_[3]) *
185  o[2];
186  break;
187  default:
188  IMP_THROW("Out of range coordinate " << coord, IndexException);
189  }
190  }
191 #endif
192  //! Rotate a vector around the origin
193  Vector3D get_rotated(const Vector3D &o) const {
194  if (!has_cache_) fill_cache();
195  return Vector3D(o * matrix_[0], o * matrix_[1], o * matrix_[2]);
196  }
197 
198 #ifndef SWIG
199  //! Get adjoint of inputs to `get_rotated` from adjoint of output
200  /** Compute the adjoint (reverse-mode sensitivity) of input vector
201  to `get_rotated` and this rotation from the adjoint of the
202  output vector.
203  */
204  void get_rotated_adjoint(const Vector3D &v, const Vector3D &Dw,
205  Vector3D *Dv, Rotation3DAdjoint *Dr) const;
206 #endif
207 
208  //! Get adjoint of inputs to `get_rotated` from adjoint of output
209  /** Compute the adjoint (reverse-mode sensitivity) of input vector
210  to `get_rotated` and this rotation from the adjoint of the
211  output vector.
212  */
213  RotatedVector3DAdjoint
214  get_rotated_adjoint(const Vector3D &v, const Vector3D &Dw) const;
215 
216  //! Get only the requested rotation coordinate of the vector
218  unsigned int coord) const {
219  if (!has_cache_) fill_cache();
220  return o * matrix_[coord];
221  }
222 
223  //! Rotate a vector around the origin
224  Vector3D operator*(const Vector3D &v) const { return get_rotated(v); }
225  Vector3D get_rotation_matrix_row(int i) const {
226  IMP_USAGE_CHECK((i >= 0) && (i <= 2), "row index out of range");
227  if (!has_cache_) fill_cache();
228  return matrix_[i];
229  }
230  IMP_SHOWABLE_INLINE(Rotation3D,
231  { out << v_[0] << " " << v_[1] << " " << v_[2] << " " << v_[3]; });
232 
233  //! Return the rotation which undoes this rotation.
234  inline Rotation3D get_inverse() const {
235  IMP_USAGE_CHECK(get_is_valid(),
236  "Attempting to invert uninitialized rotation");
237  Rotation3D ret(v_[0], -v_[1], -v_[2], -v_[3]);
238  return ret;
239  }
240 
241  //! Return the quaternion so that it can be stored
242  /** Note that there is no guarantee on which of the two
243  equivalent quaternions is returned.
244  */
245  const Vector4D &get_quaternion() const {
246  IMP_USAGE_CHECK(get_is_valid(),
247  "Attempting to access uninitialized rotation");
248  return v_;
249  }
250 
251  //! Multiply two rotations
252  Rotation3D operator*(const Rotation3D &q) const {
253  IMP_USAGE_CHECK(get_is_valid(),
254  "Attempting to compose uninitialized rotation");
255  return compose(*this, q);
256  }
257 
258  //! Compute the rotation which when composed with r gives this
259  Rotation3D operator/(const Rotation3D &r) const {
260  IMP_USAGE_CHECK(get_is_valid(),
261  "Attempting to compose uninitialized rotation");
262  return compose(*this, r.get_inverse());
263  }
264 
265  const Rotation3D &operator/=(const Rotation3D &r) {
266  *this = *this / r;
267  return *this;
268  }
269 
270  //! Return the gradient of rotated vector wrt the ith quaternion element.
271  /** Given the rotation \f$x = R(q) v\f$, where \f$v\f$ is a vector,
272  \f$q=(q_0,q_1,q_2,q_3)\f$ is the quaternion of the rotation with
273  elements \f$q_i\f$, and \f$R(q)\f$ is the corresponding rotation matrix,
274  the function returns the gradient \f$\nabla_{q_i} x\f$.
275 
276  This function just returns a single column from `get_jacobian()`, so it is
277  more efficient to call that function if all columns are needed.
278 
279  \param[in] v vector to be rotated by rotation \f$R(q)\f$
280  \param[in] i the quaternion element
281  \param[in] wrt_unnorm Gradient is computed wrt unnormalized quaternion.
282  Rotation includes a normalization operation, and
283  the gradient is projected to the tangent space at
284  \f$q\f$.
285  */
286  Vector3D get_gradient_of_rotated(const Vector3D &v, unsigned int i,
287  bool wrt_unnorm = false) const;
288 
289  //! Return the Jacobian of rotated vector wrt the quaternion.
290  /** Given the rotation \f$x = R(q) v\f$, where \f$v\f$ is a vector,
291  \f$q\f$ is the quaternion of the rotation, and \f$R(q)\f$ is the
292  corresponding rotation matrix, the function returns the 3x4 matrix
293  \f$J\f$ with elements \f$J_{ij}=\frac{\partial x_i}{\partial q_j}\f$.
294 
295  \param[in] v vector to be rotated by rotation \f$R(q)\f$
296  \param[in] wrt_unnorm Jacobian is computed wrt unnormalized quaternion.
297  Rotation includes a normalization operation, and
298  the columns are projected to the tangent space at
299  \f$q\f$.
300  */
301  Eigen::MatrixXd get_jacobian_of_rotated(
302  const Eigen::Vector3d &v, bool wrt_unnorm = false) const;
303 
304  /** Return true is the rotation is valid, false if
305  invalid or not initialized (e.g., only initialized by
306  the empty constructor)
307  */
308  bool get_is_valid() const {
309  return v_.get_squared_magnitude() > 0; // TODO: add that magnitude ~ 1?
310  }
311 };
312 
314 
315 
316 //! Get gradient of quaternion product with respect to first quaternion.
317 /** Given the rotation \f$R(p)\f$ followed by \f$R(q)\f$, where \f$p\f$ and
318  \f$q\f$ are quaternions, the quaternion of the composed rotation
319  \f$R(s)=R(q) R(p)\f$ can be expressed through the Hamilton product of the
320  two quaternions \f$s(q,p) = q p\f$. This function returns the matrix
321  \f$J\f$ with elements \f$J_{ij}=\frac{\partial s_i}{\partial q_j}\f$.
322 
323  \param[in] q rotation corresponding to first quaternion
324  \param[in] p rotation corresponding to second quaternion
325  \param[in] wrt_unnorm Jacobian is computed wrt unnormalized quaternion.
326  Rotation includes a normalization operation, and
327  the columns are projected to the tangent space at
328  \f$q\f$.
329  */
330 IMPALGEBRAEXPORT Eigen::MatrixXd get_jacobian_of_composed_wrt_first(
331  const Rotation3D &q, const Rotation3D &p, bool wrt_unnorm = false);
332 
333 //! Get Jacobian of quaternion product with respect to second quaternion.
334 /** Given the rotation \f$R(p)\f$ followed by \f$R(q)\f$, where \f$p\f$ and
335  \f$q\f$ are quaternions, the quaternion of the composed rotation
336  \f$R(s)=R(q) R(p)\f$ can be expressed through the Hamilton product of the
337  two quaternions \f$s(q,p) = q p\f$. This function returns the matrix
338  \f$J\f$ with elements \f$J_{ij}=\frac{\partial s_i}{\partial p_j}\f$.
339 
340  \param[in] q rotation corresponding to first quaternion
341  \param[in] p rotation corresponding to second quaternion
342  \param[in] wrt_unnorm Jacobian is computed wrt unnormalized quaternion.
343  Rotation includes a normalization operation, and
344  the columns are projected to the tangent space at
345  \f$p\f$.
346  */
347 IMPALGEBRAEXPORT Eigen::MatrixXd get_jacobian_of_composed_wrt_second(
348  const Rotation3D &q, const Rotation3D &p, bool wrt_unnorm = false);
349 
350 //! Return a rotation that does not do anything
351 /** \see Rotation3D */
352 inline Rotation3D get_identity_rotation_3d() { return Rotation3D(1, 0, 0, 0); }
353 
354 //! Return a distance between the two rotations
355 /** The distance runs between 0 and 1. More precisely,
356  the distance returned is distance between the two
357  quaternion vectors properly normalized, divided
358  by sqrt(2).
359 
360  A vector with distance d from the unit vector
361  represents a rotation of
362 
363  \f$ \theta= \cos^{-1}\left(1-\sqrt{2}d\right)\f$
364  \see Rotation3D
365 */
366 inline double get_distance(const Rotation3D &r0, const Rotation3D &r1) {
367  double dot =
368  (r0.get_quaternion() - r1.get_quaternion()).get_squared_magnitude();
369  double odot =
370  (r0.get_quaternion() + r1.get_quaternion()).get_squared_magnitude();
371  double ans = std::min(dot, odot);
372  // TODO: barak - added static for efficiency
373  static const double s2 = std::sqrt(2.0);
374  double ret = ans / s2;
375  return std::max(std::min(ret, 1.0), 0.0);
376 }
377 
378 //! Generate a Rotation3D object from a rotation around an axis
379 //! that is assumed to be normalized
380 /**
381  \param[in] axis_norm the normalized rotation axis passing through (0,0,0)
382  \param[in] angle the rotation angle in radians in the
383  clockwise direction
384  \note https://en.wikipedia.org/wiki/Rotation_matrix
385  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
386  angleToQuaternion/index.htm
387  \see Rotation3D
388 */
389 inline Rotation3D
391 (const Vector3D &axis_norm, double angle)
392 {
393  IMP_USAGE_CHECK(axis_norm.get_magnitude() - 1.0 < 1e-6,
394  "expected normalized vector as axis of rotation");
395  double s = std::sin(angle / 2);
396  double a, b, c, d;
397  a = std::cos(angle / 2);
398  b = axis_norm[0] * s;
399  c = axis_norm[1] * s;
400  d = axis_norm[2] * s;
401  return Rotation3D(a, b, c, d);
402 }
403 
404 
405 //! Generate a Rotation3D object from a rotation around an axis
406 /**
407  \param[in] axis the rotation axis passes through (0,0,0)
408  \param[in] angle the rotation angle in radians in the
409  clockwise direction
410  \note https://en.wikipedia.org/wiki/Rotation_matrix
411  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
412  angleToQuaternion/index.htm
413  \see Rotation3D
414 */
415 inline Rotation3D
417 (const Vector3D &axis, double angle)
418 {
419  // normalize the vector
420  Vector3D axis_norm = axis.get_unit_vector();
421  return get_rotation_about_normalized_axis(axis_norm, angle);
422 }
423 
424 //! Create a rotation from the first vector to the second one.
425 /** \see Rotation3D
426  */
427 IMPALGEBRAEXPORT Rotation3D
429 
430 //! Generate a Rotation3D object from a rotation matrix
431 /**
432  \see Rotation3D
433 */
434 IMPALGEBRAEXPORT Rotation3D
435  get_rotation_from_matrix(double m00, double m01, double m02, double m10,
436  double m11, double m12, double m20, double m21,
437  double m22);
438 
439 //! Generate a Rotation3D object from a rotation matrix
440 /**
441  \see Rotation3D
442 */
443 IMPALGEBRAEXPORT Rotation3D get_rotation_from_matrix(Eigen::Matrix3d m);
444 
445 //! Pick a rotation at random from all possible rotations
446 /** \see Rotation3D */
447 IMPALGEBRAEXPORT Rotation3D get_random_rotation_3d();
448 
449 //! Pick a rotation at random near the provided one
450 /** This method generates a rotation that is within the provided
451  distance of center.
452  \param[in] center The center of the rotational volume
453  \param[in] distance See
454  get_distance(const Rotation3D&,const Rotation3D&)
455  for a full definition.
456 
457  \note The cost of this operation increases as distance goes to 0.
458 
459  \see Rotation3D
460 */
461 IMPALGEBRAEXPORT Rotation3D
462  get_random_rotation_3d(const Rotation3D &center, double distance);
463 
464 //! Cover the space of rotations evenly
465 /** If you care about the distance between samples instead of the number
466  of samples, the "surface area" of the set of rotations is pi^2. If
467  you allocate each sample a volume of 4/3 pi d^3 (to space them d apart),
468  Then you want 3/4 pi/d^3 points.
469 
470  Creates at least num_points rotations.
471 */
472 IMPALGEBRAEXPORT Rotation3Ds
473  get_uniform_cover_rotations_3d(unsigned int num_points);
474 
475 //! Generates a nondegenerate set of Euler angles with a delta resolution
476 /**
477 \param[in] delta sample every delta angles in radians.
478  */
480  double delta);
481 
482 //! Compute a rotation from an unnormalized quaternion
483 /** \see Rotation3D */
485  VectorD<4> uv = v.get_unit_vector();
486  return Rotation3D(uv[0], uv[1], uv[2], uv[3]);
487 }
488 
489 /** \see Rotation3D
490  */
491 inline Rotation3D compose(const Rotation3D &a, const Rotation3D &b) {
492  return Rotation3D(a.v_[0] * b.v_[0] - a.v_[1] * b.v_[1] - a.v_[2] * b.v_[2] -
493  a.v_[3] * b.v_[3],
494  a.v_[0] * b.v_[1] + a.v_[1] * b.v_[0] + a.v_[2] * b.v_[3] -
495  a.v_[3] * b.v_[2],
496  a.v_[0] * b.v_[2] - a.v_[1] * b.v_[3] + a.v_[2] * b.v_[0] +
497  a.v_[3] * b.v_[1],
498  a.v_[0] * b.v_[3] + a.v_[1] * b.v_[2] - a.v_[2] * b.v_[1] +
499  a.v_[3] * b.v_[0]);
500 }
501 
502 #ifndef SWIG
503 //! Get adjoint of inputs to `compose` from adjoint of output
504 /** Compute the adjoint (reverse-mode sensitivity) of input rotations
505  to `compose` from the adjoint of the output rotation.
506  */
507 IMPALGEBRAEXPORT void
508 compose_adjoint(const Rotation3D &A, const Rotation3D &B, Vector4D DC,
509  Rotation3DAdjoint *DA, Rotation3DAdjoint *DB);
510 #endif
511 
512 //! Get adjoint of inputs to `compose` from adjoint of output
513 /** Compute the adjoint (reverse-mode sensitivity) of input rotations
514  to `compose` from the adjoint of the output rotation.
515  */
516 IMPALGEBRAEXPORT ComposeRotation3DAdjoint
517 compose_adjoint(const Rotation3D &A, const Rotation3D &B, const Rotation3DAdjoint &DC);
518 
519 /** \name Euler Angles
520  There are many conventions for how to define Euler angles, based on choices
521  of which of the x,y,z axis to use in what order and whether the rotation
522  axis is in the body frame (and hence affected by previous rotations) or in
523  in a fixed frame. See
524  https://en.wikipedia.org/wiki/Conversion_between_quaternions_and_Euler_angles
525  for a general description.
526 
527  - All Euler angles are specified in radians.
528  - Only x-y-z order is currently supported.
529  - To convert Euler angles in a different order, one can compose a Rotation3D
530  from three rotations using get_rotation_about_axis function.
531  @{
532 */
533 
534 //! Initialize a rotation in x-y-z order from three angles
535 /** \param[in] xr Rotation around the X axis in radians
536  \param[in] yr Rotation around the Y axis in radians
537  \param[in] zr Rotation around the Z axis in radians
538  \note The three rotations are represented in the original (fixed)
539  coordinate frame.
540  \see Rotation3D
541  \see FixedXYZ
542 */
543 IMPALGEBRAEXPORT Rotation3D
544  get_rotation_from_fixed_xyz(double xr, double yr, double zr);
545 
546 //! Initialize a rotation from Euler angles
547 /**
548  \param[in] phi Rotation around the Z axis in radians
549  \param[in] theta Rotation around the X axis in radians
550  \param[in] psi Rotation around the Z axis in radians
551  \note The first rotation is by an angle phi about the z-axis.
552  The second rotation is by an angle theta in [0,pi] about the
553  former x-axis , and the third rotation is by an angle psi
554  about the former z-axis.
555  \see Rotation3D
556 */
557 IMPALGEBRAEXPORT Rotation3D
558  get_rotation_from_fixed_zxz(double phi, double theta, double psi);
559 
560 //! Generate a rotation object from Euler Angles
561 /** \note The first rotation is by an angle about the z-axis.
562  The second rotation is by an angle about the new y-axis.
563  The third rotation is by an angle about the new z-axis.
564  \param[in] Rot First Euler angle (radians) defining the rotation (Z axis)
565  \param[in] Tilt Second Euler angle (radians) defining the rotation (Y axis)
566  \param[in] Psi Third Euler angle (radians) defining the rotation (Z axis)
567  \see Rotation3D
568 */
569 IMPALGEBRAEXPORT Rotation3D
570  get_rotation_from_fixed_zyz(double Rot, double Tilt, double Psi);
571 
572 //! A simple class for returning XYZ Euler angles
573 class FixedXYZ : public GeometricPrimitiveD<3> {
574  double v_[3];
575 
576  friend class cereal::access;
577 
578  template<class Archive> void serialize(Archive &ar) {
579  ar(v_[0], v_[1], v_[2]);
580  }
581 
582  public:
583  FixedXYZ() {}
584  FixedXYZ(double x, double y, double z) {
585  v_[0] = x;
586  v_[1] = y;
587  v_[2] = z;
588  }
589  double get_x() const { return v_[0]; }
590  double get_y() const { return v_[1]; }
591  double get_z() const { return v_[2]; }
593  { out << v_[0] << " " << v_[1] << " " << v_[2]; });
594 };
595 
597 
598 //! The inverse of rotation_from_fixed_xyz()
599 /**
600  \see rotation_from_fixed_xyz()
601  \see Rotation3D
602  \see FixesXYZ
603 */
604 IMPALGEBRAEXPORT FixedXYZ get_fixed_xyz_from_rotation(const Rotation3D &r);
605 
606 /** @}*/
607 
608 //! Interpolate between two rotations
609 /** It f ==0, return b, if f==1 return a.
610  \see Rotation3D */
612  double f) {
613  VectorD<4> bq = b.get_quaternion(), aq = a.get_quaternion();
614  if (bq * aq < 0) bq = -bq;
615  return Rotation3D(f * aq + (1 - f) * bq);
616 }
617 
618 /** Return the rotation which takes the native x and y axes to the
619  given x and y axes.
620  The two axes must be perpendicular unit vectors.
621 */
622 IMPALGEBRAEXPORT Rotation3D
623  get_rotation_from_x_y_axes(const Vector3D &x, const Vector3D &y);
624 
625 //! Decompose a Rotation3D object into a rotation around an axis
626 /** For all identity rotations, returns the axis [1,0,0] and the angle 0.0.
627 
628  \note https://en.wikipedia.org/wiki/Rotation_matrix
629  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
630  angleToQuaternion/index.htm
631  \see Rotation3D
632 
633  @return axis direction and rotation about the axis in Radians
634 */
635 IMPALGEBRAEXPORT std::pair<Vector3D, double> get_axis_and_angle(
636  const Rotation3D &rot);
637 
638 typedef std::pair<Vector3D, double> AxisAnglePair;
639 #ifndef IMP_DOXYGEN
640 typedef Vector<AxisAnglePair> AxisAnglePairs;
641 #endif
642 
643 IMPALGEBRA_END_NAMESPACE
644 #endif /* IMPALGEBRA_ROTATION_3D_H */
Vector3D get_rotated(const Vector3D &o) const
Rotate a vector around the origin.
Definition: Rotation3D.h:193
Rotation3D get_rotation_about_normalized_axis(const Vector3D &axis_norm, double angle)
Definition: Rotation3D.h:391
Rotation3D(const Rotation3D &rot)
Rotation3D copy constructor.
Definition: Rotation3D.h:106
Base class for geometric types.
#define IMP_SHOWABLE_INLINE(Name, how_to_show)
Declare the methods needed by an object that can be printed.
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:125
#define IMP_USAGE_CHECK_FLOAT_EQUAL(expra, exprb, message)
Definition: check_macros.h:178
Rotation2D compose(const Rotation2D &a, const Rotation2D &b)
Compose two rotations a and b.
Definition: Rotation2D.h:121
algebra::Rotation3Ds get_uniformly_sampled_rotations(double delta)
Generates a nondegenerate set of Euler angles with a delta resolution.
Rotation3D get_random_rotation_3d(const Rotation3D &center, double distance)
Pick a rotation at random near the provided one.
Rotation3D operator/(const Rotation3D &r) const
Compute the rotation which when composed with r gives this.
Definition: Rotation3D.h:259
A simple class for returning XYZ Euler angles.
Definition: Rotation3D.h:573
Rotation3D(const VectorD< 4 > &v, bool assume_normalized=false)
Definition: Rotation3D.h:119
Vector3D operator*(const Vector3D &v) const
Rotate a vector around the origin.
Definition: Rotation3D.h:224
Rotation3D compose(const Rotation3D &a, const Rotation3D &b)
Definition: Rotation3D.h:491
double get_rotated_one_coordinate(const Vector3D &o, unsigned int coord) const
Get only the requested rotation coordinate of the vector.
Definition: Rotation3D.h:217
VT get_unit_vector(VT vt)
Returns a unit vector pointing at the same direction as this vector.
Definition: VectorBaseD.h:272
Eigen::MatrixXd get_jacobian_of_composed_wrt_first(const Rotation3D &q, const Rotation3D &p, bool wrt_unnorm=false)
Get gradient of quaternion product with respect to first quaternion.
bool get_is_valid() const
Definition: Rotation3D.h:308
const Vector4D & get_quaternion() const
Return the quaternion so that it can be stored.
Definition: Rotation3D.h:245
VectorD< 4 > Vector4D
Definition: VectorD.h:412
Rotation3D get_rotation_about_axis(const Vector3D &axis, double angle)
Generate a Rotation3D object from a rotation around an axis.
Definition: Rotation3D.h:417
#define IMP_VALUES(Name, PluralName)
Define the type for storing sets of values.
Definition: value_macros.h:23
Base class for geometric types.
Functions to deal with very common math operations.
A Cartesian vector in D-dimensions.
Definition: VectorD.h:39
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:234
Rotation3D get_rotation_from_fixed_zyz(double Rot, double Tilt, double Psi)
Generate a rotation object from Euler Angles.
std::pair< Vector3D, double > get_axis_and_angle(const Rotation3D &rot)
Decompose a Rotation3D object into a rotation around an axis.
double get_distance(const Rotation3D &r0, const Rotation3D &r1)
Return a distance between the two rotations.
Definition: Rotation3D.h:366
Rotation3D get_interpolated(const Rotation3D &a, const Rotation3D &b, double f)
Interpolate between two rotations.
Definition: Rotation3D.h:611
#define IMP_NO_SWIG(x)
Hide the line when SWIG is compiled or parses it.
Definition: swig_macros.h:18
Rotation3D get_rotation_from_vector4d(const VectorD< 4 > &v)
Compute a rotation from an unnormalized quaternion.
Definition: Rotation3D.h:484
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:252
3D rotation class.
Definition: Rotation3D.h:52
Eigen::MatrixXd get_jacobian_of_composed_wrt_second(const Rotation3D &q, const Rotation3D &p, bool wrt_unnorm=false)
Get Jacobian of quaternion product with respect to second quaternion.
#define IMP_THROW(message, exception_name)
Throw an exception with a message.
Definition: check_macros.h:50
Rotation3D(double a, double b, double c, double d)
Create a rotation from a quaternion.
Definition: Rotation3D.h:129
VectorD< 3 > Vector3D
Definition: VectorD.h:408
Rotation3D get_identity_rotation_3d()
Return a rotation that does not do anything.
Definition: Rotation3D.h:352
ComposeRotation3DAdjoint compose_adjoint(const Rotation3D &A, const Rotation3D &B, const Rotation3DAdjoint &DC)
Get adjoint of inputs to compose from adjoint of output.
Simple 3D vector class.
#define IMP_USAGE_CHECK(expr, message)
A runtime test for incorrect usage of a class or method.
Definition: check_macros.h:168
Logging and error reporting support.
Rotation3D get_rotation_from_fixed_zxz(double phi, double theta, double psi)
Initialize a rotation from Euler angles.
Rotation3D get_rotation_from_matrix(Eigen::Matrix3d m)
Generate a Rotation3D object from a rotation matrix.
IMP::Vector< Rotation3D > Rotation3Ds
Definition: Rotation3D.h:313