8 #ifndef IMPALGEBRA_ROTATION_3D_H
9 #define IMPALGEBRA_ROTATION_3D_H
11 #include <IMP/algebra/algebra_config.h>
16 #include <Eigen/Dense>
23 IMPALGEBRA_BEGIN_NAMESPACE
25 #if !defined(IMP_DOXYGEN) && !defined(SWIG)
27 Rotation3D
compose(
const Rotation3D &a,
const Rotation3D &b);
48 mutable bool has_cache_;
52 void fill_cache()
const {
54 "Attempting to apply uninitialized rotation");
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];
67 Vector3D(v0s + v1s - v2s - v3s, 2 * (v12 - v03), 2 * (v13 + v02));
69 Vector3D(2 * (v12 + v03), v0s - v1s + v2s - v3s, 2 * (v23 - v01));
71 Vector3D(2 * (v13 - v02), 2 * (v23 + v01), v0s - v1s - v2s + v3s);
78 v_(rot.v_), has_cache_(rot.has_cache_)
81 matrix_[0]=rot.matrix_[0];
82 matrix_[1]=rot.matrix_[1];
83 matrix_[2]=rot.matrix_[2];
93 bool assume_normalized=
false)
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());
120 "Attempting to access uninitialized rotation");
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]) *
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]) *
136 double get_rotated_one_coordinate_no_cache(
const Vector3D &o,
137 unsigned int coord)
const {
139 "Attempting to apply uninitialized rotation");
142 return (v_[0] * v_[0] + v_[1] * v_[1] - v_[2] * v_[2] - v_[3] * v_[3]) *
144 2 * (v_[1] * v_[2] - v_[0] * v_[3]) * o[1] +
145 2 * (v_[1] * v_[3] + v_[0] * v_[2]) * o[2];
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]) *
151 2 * (v_[2] * v_[3] - v_[0] * v_[1]) * o[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]) *
161 IMP_THROW(
"Out of range coordinate " << coord, IndexException);
167 if (!has_cache_) fill_cache();
168 return Vector3D(o * matrix_[0], o * matrix_[1], o * matrix_[2]);
173 unsigned int coord)
const {
174 if (!has_cache_) fill_cache();
175 return o * matrix_[coord];
180 Vector3D get_rotation_matrix_row(
int i)
const {
182 if (!has_cache_) fill_cache();
186 { out << v_[0] <<
" " << v_[1] <<
" " << v_[2] <<
" " << v_[3]; });
191 "Attempting to invert uninitialized rotation");
192 Rotation3D ret(v_[0], -v_[1], -v_[2], -v_[3]);
202 "Attempting to access uninitialized rotation");
209 "Attempting to compose uninitialized rotation");
216 "Attempting to compose uninitialized rotation");
240 bool projected =
true)
const;
253 Eigen::MatrixXd get_gradient(
254 const Eigen::Vector3d &v,
bool projected =
true)
const;
261 return v_.get_squared_magnitude() > 0;
321 double ans = std::min(dot, odot);
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);
344 "expected normalized vector as axis of rotation");
345 double s = std::sin(angle / 2);
347 a = std::cos(angle / 2);
348 b = axis_norm[0] * s;
349 c = axis_norm[1] * s;
350 d = axis_norm[2] * s;
370 Vector3D axis_norm = axis.get_unit_vector();
377 IMPALGEBRAEXPORT Rotation3D
384 IMPALGEBRAEXPORT Rotation3D
386 double m11,
double m12,
double m20,
double m21,
411 IMPALGEBRAEXPORT Rotation3D
436 return Rotation3D(uv[0], uv[1], uv[2], uv[3]);
442 return Rotation3D(a.v_[0] * b.v_[0] - a.v_[1] * b.v_[1] - a.v_[2] * b.v_[2] -
444 a.v_[0] * b.v_[1] + a.v_[1] * b.v_[0] + a.v_[2] * b.v_[3] -
446 a.v_[0] * b.v_[2] - a.v_[1] * b.v_[3] + a.v_[2] * b.v_[0] +
448 a.v_[0] * b.v_[3] + a.v_[1] * b.v_[2] - a.v_[2] * b.v_[1] +
476 IMPALGEBRAEXPORT Rotation3D
490 IMPALGEBRAEXPORT Rotation3D
502 IMPALGEBRAEXPORT Rotation3D
511 FixedXYZ(
double x,
double y,
double z) {
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]; });
541 if (bq * aq < 0) bq = -bq;
549 IMPALGEBRAEXPORT Rotation3D
563 const Rotation3D &rot);
565 typedef std::pair<Vector3D, double> AxisAnglePair;
570 IMPALGEBRA_END_NAMESPACE
Vector3D get_rotated(const Vector3D &o) const
Rotate a vector around the origin.
Rotation3D get_rotation_about_normalized_axis(const Vector3D &axis_norm, double angle)
Rotation3D(const Rotation3D &rot)
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.
#define IMP_USAGE_CHECK_FLOAT_EQUAL(expra, exprb, message)
Rotation2D compose(const Rotation2D &a, const Rotation2D &b)
Compose two rotations a and b.
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 ¢er, 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.
A simple class for returning XYZ Euler angles.
Rotation3D(const VectorD< 4 > &v, bool assume_normalized=false)
IMP_CXX11_DEFAULT_COPY_CONSTRUCTOR(Rotation3D);.
Vector3D operator*(const Vector3D &v) const
Rotate a vector around the origin.
Rotation3D compose(const Rotation3D &a, const Rotation3D &b)
double get_rotated_one_coordinate(const Vector3D &o, unsigned int coord) const
Get only the requested rotation coordinate of the vector.
VT get_unit_vector(VT vt)
Returns a unit vector pointing at the same direction as this vector.
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
const Vector4D & get_quaternion() const
Return the quaternion so that it can be stored.
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.
Base class for geometric types.
Functions to deal with very common math operations.
A Cartesian vector in D-dimensions.
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.
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.
Rotation3D get_interpolated(const Rotation3D &a, const Rotation3D &b, double f)
Interpolate between two rotations.
#define IMP_NO_SWIG(x)
Hide the line when SWIG is compiled or parses it.
Rotation3D get_rotation_from_vector4d(const VectorD< 4 > &v)
Compute a rotation from an unnormalized quaternion.
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.
#define IMP_THROW(message, exception_name)
Throw an exception with a message.
Rotation3D(double a, double b, double c, double d)
Create a rotation from a quaternion.
Rotation3D get_identity_rotation_3d()
Return a rotation that does not do anything.
#define IMP_USAGE_CHECK(expr, message)
A runtime test for incorrect usage of a class or method.
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