IMP logo
IMP Reference Guide  2.5.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-2015 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 <IMP/algebra/eigen3/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  IMP_CXX11_DEFAULT_COPY_CONSTRUCTOR(Rotation3D);
76 
77  //! Create a rotation from a vector of 4 quaternion coefficients.
78  //! @note: use assume_normalize with care - inputting an unnormalized
79  //! vector would result in unexpected results if it is true
80  explicit Rotation3D(const VectorD<4> &v,
81  bool assume_normalized=false)
82  : v_(assume_normalized ? v : v.get_unit_vector()),
83  has_cache_(false) {}
84 
85  //! Create an invalid rotation
86  Rotation3D() : v_(0, 0, 0, 0), has_cache_(false) {}
87  //! Create a rotation from a quaternion
88  /** \throw ValueException if the rotation is not a unit vector.
89  */
90  Rotation3D(double a, double b, double c, double d)
91  : v_(a, b, c, d), has_cache_(false) {
93  v_.get_squared_magnitude(), 1.0,
94  "Attempting to construct a rotation from a "
95  << " non-quaternion value. The coefficient vector"
96  << " must have a length of 1. Got: " << a << " " << b << " " << c
97  << " " << d << " gives " << v_.get_squared_magnitude());
98  if (a < 0) {
99  // make them canonical
100  v_ = -v_;
101  }
102  }
103  ~Rotation3D(){}
104 
105 #ifndef IMP_DOXYGEN
106  Vector3D get_rotated_no_cache(const Vector3D &o) const {
107  IMP_USAGE_CHECK(get_is_valid(),
108  "Attempting to access uninitialized rotation");
109  return Vector3D(
110  (v_[0] * v_[0] + v_[1] * v_[1] - v_[2] * v_[2] - v_[3] * v_[3]) * o[0] +
111  2 * (v_[1] * v_[2] - v_[0] * v_[3]) * o[1] +
112  2 * (v_[1] * v_[3] + v_[0] * v_[2]) * o[2],
113  2 * (v_[1] * v_[2] + v_[0] * v_[3]) * o[0] +
114  (v_[0] * v_[0] - v_[1] * v_[1] + v_[2] * v_[2] - v_[3] * v_[3]) *
115  o[1] +
116  2 * (v_[2] * v_[3] - v_[0] * v_[1]) * o[2],
117  2 * (v_[1] * v_[3] - v_[0] * v_[2]) * o[0] +
118  2 * (v_[2] * v_[3] + v_[0] * v_[1]) * o[1] +
119  (v_[0] * v_[0] - v_[1] * v_[1] - v_[2] * v_[2] + v_[3] * v_[3]) *
120  o[2]);
121  }
122 
123  //! Gets only the requested rotation coordinate of the vector
124  double get_rotated_one_coordinate_no_cache(const Vector3D &o,
125  unsigned int coord) const {
126  IMP_USAGE_CHECK(get_is_valid(),
127  "Attempting to apply uninitialized rotation");
128  switch (coord) {
129  case 0:
130  return (v_[0] * v_[0] + v_[1] * v_[1] - v_[2] * v_[2] - v_[3] * v_[3]) *
131  o[0] +
132  2 * (v_[1] * v_[2] - v_[0] * v_[3]) * o[1] +
133  2 * (v_[1] * v_[3] + v_[0] * v_[2]) * o[2];
134  break;
135  case 1:
136  return 2 * (v_[1] * v_[2] + v_[0] * v_[3]) * o[0] +
137  (v_[0] * v_[0] - v_[1] * v_[1] + v_[2] * v_[2] - v_[3] * v_[3]) *
138  o[1] +
139  2 * (v_[2] * v_[3] - v_[0] * v_[1]) * o[2];
140 
141  break;
142  case 2:
143  return 2 * (v_[1] * v_[3] - v_[0] * v_[2]) * o[0] +
144  2 * (v_[2] * v_[3] + v_[0] * v_[1]) * o[1] +
145  (v_[0] * v_[0] - v_[1] * v_[1] - v_[2] * v_[2] + v_[3] * v_[3]) *
146  o[2];
147  break;
148  default:
149  IMP_THROW("Out of range coordinate " << coord, IndexException);
150  }
151  }
152 #endif
153  //! Rotate a vector around the origin
154  Vector3D get_rotated(const Vector3D &o) const {
155  IMP_USAGE_CHECK(get_is_valid(),
156  "Attempting to apply uninitialized rotation");
157  if (!has_cache_) fill_cache();
158  return Vector3D(o * matrix_[0], o * matrix_[1], o * matrix_[2]);
159  }
160 
161  //! Gets only the requested rotation coordinate of the vector
163  unsigned int coord) const {
164  IMP_USAGE_CHECK(get_is_valid(),
165  "Attempting to apply uninitialized rotation");
166  if (!has_cache_) fill_cache();
167  return o * matrix_[coord];
168  }
169 
170  //! Rotate a vector around the origin
171  Vector3D operator*(const Vector3D &v) const { return get_rotated(v); }
172  Vector3D get_rotation_matrix_row(int i) const {
173  IMP_USAGE_CHECK((i >= 0) && (i <= 2), "row index out of range");
174  if (!has_cache_) fill_cache();
175  return matrix_[i];
176  }
177  IMP_SHOWABLE_INLINE(Rotation3D,
178  { out << v_[0] << " " << v_[1] << " " << v_[2] << " " << v_[3]; });
179 
180  //! Return the rotation which undoes this rotation.
181  inline Rotation3D get_inverse() const {
182  IMP_USAGE_CHECK(get_is_valid(),
183  "Attempting to invert uninitialized rotation");
184  Rotation3D ret(v_[0], -v_[1], -v_[2], -v_[3]);
185  return ret;
186  }
187 
188  //! Return the quaternion so that it can be stored
189  /** Note that there is no guarantee on which of the two
190  equivalent quaternions is returned.
191  */
192  const Vector4D &get_quaternion() const {
193  IMP_USAGE_CHECK(get_is_valid(),
194  "Attempting to access uninitialized rotation");
195  return v_;
196  }
197 
198  //! multiply two rotations
199  Rotation3D operator*(const Rotation3D &q) const {
200  IMP_USAGE_CHECK(get_is_valid(),
201  "Attempting to compose uninitialized rotation");
202  return compose(*this, q);
203  }
204 
205  //! Compute the rotation which when composed with r gives this
206  Rotation3D operator/(const Rotation3D &r) const {
207  IMP_USAGE_CHECK(get_is_valid(),
208  "Attempting to compose uninitialized rotation");
209  return compose(*this, r.get_inverse());
210  }
211 
212  const Rotation3D &operator/=(const Rotation3D &r) {
213  *this = *this / r;
214  return *this;
215  }
216 
217  /** \brief Return the derivative of the local position o with respect to
218  the i'th internal quaternion coefficient, for i in [0..3],
219  namely (dx/dQi, dy/dQi, dz/dQi)
220  */
221  const Vector3D get_derivative(const Vector3D &o, unsigned int i) const;
222 
223  /** returns true is the rotation is valid, false if
224  invalid or not initialized (e.g., only initialized by
225  the empty constructor)
226  */
227  bool get_is_valid() const {
228  return v_.get_squared_magnitude() > 0; // TODO: add that magnitude ~ 1?
229  }
230 };
231 
233 
234 //! Return a rotation that does not do anything
235 /** \see Rotation3D */
236 inline Rotation3D get_identity_rotation_3d() { return Rotation3D(1, 0, 0, 0); }
237 
238 //! Return a distance between the two rotations
239 /** The distance runs between 0 and 1. More precisely,
240  the distance returned is distance between the two
241  quaternion vectors properly normalized, divided
242  by sqrt(2).
243 
244  A vector with distance d from the unit vector
245  represents a rotation of
246 
247  \f$ \theta= \cos^{-1}\left(1-\sqrt{2}d\right)\f$
248  \see Rotation3D
249 */
250 inline double get_distance(const Rotation3D &r0, const Rotation3D &r1) {
251  double dot =
252  (r0.get_quaternion() - r1.get_quaternion()).get_squared_magnitude();
253  double odot =
254  (r0.get_quaternion() + r1.get_quaternion()).get_squared_magnitude();
255  double ans = std::min(dot, odot);
256  // TODO: barak - added static for efficiency
257  static const double s2 = std::sqrt(2.0);
258  double ret = ans / s2;
259  return std::max(std::min(ret, 1.0), 0.0);
260 }
261 
262 //! Generate a Rotation3D object from a rotation around an axis
263 //! that is assumed to be normalized
264 /**
265  \param[in] axis_norm the normalized rotation axis passing through (0,0,0)
266  \param[in] angle the rotation angle in radians in the
267  clockwise direction
268  \note http://en.wikipedia.org/wiki/Rotation_matrix
269  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
270  angleToQuaternion/index.htm
271  \see Rotation3D
272 */
273 IMPALGEBRAEXPORT Rotation3D
274  get_rotation_about_normalized_axis(const Vector3D &axis_norm, double angle);
275 
276 //! Generate a Rotation3D object from a rotation around an axis
277 /**
278  \param[in] axis the rotation axis passes through (0,0,0)
279  \param[in] angle the rotation angle in radians in the
280  clockwise direction
281  \note http://en.wikipedia.org/wiki/Rotation_matrix
282  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
283  angleToQuaternion/index.htm
284  \see Rotation3D
285 */
286 IMPALGEBRAEXPORT Rotation3D
287  get_rotation_about_axis(const Vector3D &axis, double angle);
288 
289 //! Create a rotation from the first vector to the second one.
290 /** \see Rotation3D
291  */
292 IMPALGEBRAEXPORT Rotation3D
294 
295 //! Generate a Rotation3D object from a rotation matrix
296 /**
297  \see Rotation3D
298 */
299 IMPALGEBRAEXPORT Rotation3D
300  get_rotation_from_matrix(double m00, double m01, double m02, double m10,
301  double m11, double m12, double m20, double m21,
302  double m22);
303 
304 //! Generate a Rotation3D object from a rotation matrix
305 /**
306  \see Rotation3D
307 */
308 IMPALGEBRAEXPORT Rotation3D get_rotation_from_matrix(IMP_Eigen::Matrix3d m);
309 
310 //! Pick a rotation at random from all possible rotations
311 /** \see Rotation3D */
312 IMPALGEBRAEXPORT Rotation3D get_random_rotation_3d();
313 
314 //! Pick a rotation at random near the provided one
315 /** This method generates a rotation that is within the provided
316  distance of center.
317  \param[in] center The center of the rotational volume
318  \param[in] distance See
319  get_distance(const Rotation3D&,const Rotation3D&)
320  for a full definition.
321 
322  \note The cost of this operation increases as distance goes to 0.
323 
324  \see Rotation3D
325 */
326 IMPALGEBRAEXPORT Rotation3D
327  get_random_rotation_3d(const Rotation3D &center, double distance);
328 
329 //! Cover the space of rotations evenly
330 /** If you care about the distance between samples instead of the number
331  of samples, the "surface area" of the set of rotations is pi^2. If
332  you allocate each sample a volume of 4/3 pi d^3 (to space them d apart),
333  Then you want 3/4 pi/d^3 points.
334 
335  Creates at least num_points rotations.
336 */
337 IMPALGEBRAEXPORT Rotation3Ds
338  get_uniform_cover_rotations_3d(unsigned int num_points);
339 
340 //! Generates a nondegenerate set of Euler angles with a delta resolution
341 /**
342 \param[in] delta sample every delta angles in radians.
343  */
345  double delta);
346 
347 //! Compute a rotation from an unnormalized quaternion
348 /** \see Rotation3D */
350  VectorD<4> uv = v.get_unit_vector();
351  return Rotation3D(uv[0], uv[1], uv[2], uv[3]);
352 }
353 
354 /** \see Rotation3D
355  */
356 inline Rotation3D compose(const Rotation3D &a, const Rotation3D &b) {
357  return Rotation3D(a.v_[0] * b.v_[0] - a.v_[1] * b.v_[1] - a.v_[2] * b.v_[2] -
358  a.v_[3] * b.v_[3],
359  a.v_[0] * b.v_[1] + a.v_[1] * b.v_[0] + a.v_[2] * b.v_[3] -
360  a.v_[3] * b.v_[2],
361  a.v_[0] * b.v_[2] - a.v_[1] * b.v_[3] + a.v_[2] * b.v_[0] +
362  a.v_[3] * b.v_[1],
363  a.v_[0] * b.v_[3] + a.v_[1] * b.v_[2] - a.v_[2] * b.v_[1] +
364  a.v_[3] * b.v_[0]);
365 }
366 
367 /** \name Euler Angles
368  There are many conventions for how to define Euler angles, based on choices
369  of which of the x,y,z axis to use in what order and whether the rotation
370  axis is in the body frame (and hence affected by previous rotations) or in
371  in a fixed frame. See
372  http://en.wikipedia.org/wiki/Conversion_between_quaternions_and_Euler_angles
373  for a general description.
374 
375  - All Euler angles are specified in radians.
376  - Only x-y-z order is currently supported.
377  - To convert Euler angles in a different order, one can compose a Rotation3D
378  from three rotations using get_rotation_about_axis function.
379  @{
380 */
381 
382 //! Initialize a rotation in x-y-z order from three angles
383 /** \param[in] xr Rotation around the X axis in radians
384  \param[in] yr Rotation around the Y axis in radians
385  \param[in] zr Rotation around the Z axis in radians
386  \note The three rotations are represented in the original (fixed)
387  coordinate frame.
388  \see Rotation3D
389  \see FixedXYZ
390 */
391 IMPALGEBRAEXPORT Rotation3D
392  get_rotation_from_fixed_xyz(double xr, double yr, double zr);
393 
394 //! Initialize a rotation from Euler angles
395 /**
396  \param[in] phi Rotation around the Z axis in radians
397  \param[in] theta Rotation around the X axis in radians
398  \param[in] psi Rotation around the Z axis in radians
399  \note The first rotation is by an angle phi about the z-axis.
400  The second rotation is by an angle theta in [0,pi] about the
401  former x-axis , and the third rotation is by an angle psi
402  about the former z-axis.
403  \see Rotation3D
404 */
405 IMPALGEBRAEXPORT Rotation3D
406  get_rotation_from_fixed_zxz(double phi, double theta, double psi);
407 
408 //! Generate a rotation object from Euler Angles
409 /** \note The first rotation is by an angle about the z-axis.
410  The second rotation is by an angle about the new y-axis.
411  The third rotation is by an angle about the new z-axis.
412  \param[in] Rot First Euler angle (radians) defining the rotation (Z axis)
413  \param[in] Tilt Second Euler angle (radians) defining the rotation (Y axis)
414  \param[in] Psi Third Euler angle (radians) defining the rotation (Z axis)
415  \see Rotation3D
416 */
417 IMPALGEBRAEXPORT Rotation3D
418  get_rotation_from_fixed_zyz(double Rot, double Tilt, double Psi);
419 
420 //! A simple class for returning XYZ Euler angles
421 class FixedXYZ : public GeometricPrimitiveD<3> {
422  double v_[3];
423 
424  public:
425  FixedXYZ() {}
426  FixedXYZ(double x, double y, double z) {
427  v_[0] = x;
428  v_[1] = y;
429  v_[2] = z;
430  }
431  double get_x() const { return v_[0]; }
432  double get_y() const { return v_[1]; }
433  double get_z() const { return v_[2]; }
435  { out << v_[0] << " " << v_[1] << " " << v_[2]; });
436 };
437 
439 
440 //! The inverse of rotation_from_fixed_xyz()
441 /**
442  \see rotation_from_fixed_xyz()
443  \see Rotation3D
444  \see FixesXYZ
445 */
446 IMPALGEBRAEXPORT FixedXYZ get_fixed_xyz_from_rotation(const Rotation3D &r);
447 
448 /** @}*/
449 
450 //! Interpolate between two rotations
451 /** It f ==0, return b, if f==1 return a.
452  \see Rotation3D */
454  double f) {
455  VectorD<4> bq = b.get_quaternion(), aq = a.get_quaternion();
456  if (bq * aq < 0) bq = -bq;
457  return Rotation3D(f * aq + (1 - f) * bq);
458 }
459 
460 /** Return the rotation which takes the native x and y axes to the
461  given x and y axes.
462  The two axes must be perpendicular unit vectors.
463 */
464 IMPALGEBRAEXPORT Rotation3D
465  get_rotation_from_x_y_axes(const Vector3D &x, const Vector3D &y);
466 
467 //! Decompose a Rotation3D object into a rotation around an axis
468 /** For all identity rotations, returns the axis [1,0,0] and the angle 0.0.
469 
470  \note http://en.wikipedia.org/wiki/Rotation_matrix
471  \note www.euclideanspace.com/maths/geometry/rotations/conversions/
472  angleToQuaternion/index.htm
473  \see Rotation3D
474 
475  @return axis direction and rotation about the axis in Radians
476 */
477 IMPALGEBRAEXPORT std::pair<Vector3D, double> get_axis_and_angle(
478  const Rotation3D &rot);
479 
480 typedef std::pair<Vector3D, double> AxisAnglePair;
481 #ifndef IMP_DOXYGEN
482 typedef Vector<AxisAnglePair> AxisAnglePairs;
483 #endif
484 
485 IMPALGEBRA_END_NAMESPACE
486 #endif /* IMPALGEBRA_ROTATION_3D_H */
Vector3D get_rotated(const Vector3D &o) const
Rotate a vector around the origin.
Definition: Rotation3D.h:154
Rotation3D get_rotation_about_normalized_axis(const Vector3D &axis_norm, double angle)
Basic types used by IMP.
#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:86
#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:206
A simple class for returning XYZ Euler angles.
Definition: Rotation3D.h:421
Rotation3D(const VectorD< 4 > &v, bool assume_normalized=false)
Definition: Rotation3D.h:80
Vector3D operator*(const Vector3D &v) const
Rotate a vector around the origin.
Definition: Rotation3D.h:171
Rotation3D compose(const Rotation3D &a, const Rotation3D &b)
Definition: Rotation3D.h:356
double get_rotated_one_coordinate(const Vector3D &o, unsigned int coord) const
Gets only the requested rotation coordinate of the vector.
Definition: Rotation3D.h:162
VT get_unit_vector(VT vt)
Returns a unit vector pointing at the same direction as this vector.
Definition: VectorBaseD.h:229
bool get_is_valid() const
Definition: Rotation3D.h:227
const Vector4D & get_quaternion() const
Return the quaternion so that it can be stored.
Definition: Rotation3D.h:192
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.
#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:181
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:250
Rotation3D get_interpolated(const Rotation3D &a, const Rotation3D &b, double f)
Interpolate between two rotations.
Definition: Rotation3D.h:453
#define IMP_NO_SWIG(x)
Hide the line when SWIG is compiled or parses it.
Definition: swig_macros.h:18
Rotation3D get_rotation_from_matrix(IMP_Eigen::Matrix3d m)
Generate a Rotation3D object from a rotation matrix.
Rotation3D get_rotation_from_vector4d(const VectorD< 4 > &v)
Compute a rotation from an unnormalized quaternion.
Definition: Rotation3D.h:349
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.
Rotation3D operator*(const Rotation3D &q) const
multiply two rotations
Definition: Rotation3D.h:199
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:90
VectorD< 3 > Vector3D
Definition: VectorD.h:395
Rotation3D get_identity_rotation_3d()
Return a rotation that does not do anything.
Definition: Rotation3D.h:236
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.
IMP::Vector< Rotation3D > Rotation3Ds
Definition: Rotation3D.h:232