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