Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
Quaternion.h
Go to the documentation of this file.
1// @(#)root/mathcore:$Id$
2// Authors: W. Brown, M. Fischler, L. Moneta 2005
3
4 /**********************************************************************
5 * *
6 * Copyright (c) 2005 , LCG ROOT FNAL MathLib Team *
7 * *
8 * *
9 **********************************************************************/
10
11// Header file for rotation in 3 dimensions, represented by a quaternion
12// Created by: Mark Fischler Thurs June 9 2005
13//
14// Last update: $Id$
15//
16#ifndef ROOT_Math_GenVector_Quaternion
17#define ROOT_Math_GenVector_Quaternion 1
18
19
26
27#include <algorithm>
28#include <cassert>
29
30
31namespace ROOT {
32namespace Math {
33
34
35//__________________________________________________________________________________________
36 /**
37 Rotation class with the (3D) rotation represented by
38 a unit quaternion (u, i, j, k).
39 This is the optimal representation for multiplication of multiple
40 rotations, and for computation of group-manifold-invariant distance
41 between two rotations.
42 See also ROOT::Math::AxisAngle, ROOT::Math::EulerAngles, and ROOT::Math::Rotation3D.
43
44 @ingroup GenVector
45 */
46
48
49public:
50
51 typedef double Scalar;
52
53 // ========== Constructors and Assignment =====================
54
55 /**
56 Default constructor (identity rotation)
57 */
59 : fU(1.0)
60 , fI(0.0)
61 , fJ(0.0)
62 , fK(0.0)
63 { }
64
65 /**
66 Construct given a pair of pointers or iterators defining the
67 beginning and end of an array of four Scalars
68 */
69 template<class IT>
70 Quaternion(IT begin, IT end) { SetComponents(begin,end); }
71
72 // ======== Construction From other Rotation Forms ==================
73
74 /**
75 Construct from another supported rotation type (see gv_detail::convert )
76 */
77 template <class OtherRotation>
78 explicit Quaternion(const OtherRotation & r) {gv_detail::convert(r,*this);}
79
80
81 /**
82 Construct from four Scalars representing the coefficients of u, i, j, k
83 */
85 fU(u), fI(i), fJ(j), fK(k) { }
86
87 // The compiler-generated copy ctor, copy assignment, and dtor are OK.
88
89 /**
90 Re-adjust components to eliminate small deviations from |Q| = 1
91 orthonormality.
92 */
93 void Rectify();
94
95 /**
96 Assign from another supported rotation type (see gv_detail::convert )
97 */
98 template <class OtherRotation>
99 Quaternion & operator=( OtherRotation const & r ) {
100 gv_detail::convert(r,*this);
101 return *this;
102 }
103
104 // ======== Components ==============
105
106 /**
107 Set the four components given an iterator to the start of
108 the desired data, and another to the end (4 past start).
109 */
110 template<class IT>
111 void SetComponents(IT begin, IT end) {
112 fU = *begin++;
113 fI = *begin++;
114 fJ = *begin++;
115 fK = *begin++;
116 (void)end;
117 assert (end==begin);
118 }
119
120 /**
121 Get the components into data specified by an iterator begin
122 and another to the end of the desired data (4 past start).
123 */
124 template<class IT>
125 void GetComponents(IT begin, IT end) const {
126 *begin++ = fU;
127 *begin++ = fI;
128 *begin++ = fJ;
129 *begin++ = fK;
130 (void)end;
131 assert (end==begin);
132 }
133
134 /**
135 Get the components into data specified by an iterator begin
136 */
137 template<class IT>
138 void GetComponents(IT begin ) const {
139 *begin++ = fU;
140 *begin++ = fI;
141 *begin++ = fJ;
142 *begin = fK;
143 }
144
145 /**
146 Set the components based on four Scalars. The sum of the squares of
147 these Scalars should be 1; no checking is done.
148 */
150 fU=u; fI=i; fJ=j; fK=k;
151 }
152
153 /**
154 Get the components into four Scalars.
155 */
156 void GetComponents(Scalar & u, Scalar & i, Scalar & j, Scalar & k) const {
157 u=fU; i=fI; j=fJ; k=fK;
158 }
159
160 /**
161 Access to the four quaternion components:
162 U() is the coefficient of the identity Pauli matrix,
163 I(), J() and K() are the coefficients of sigma_x, sigma_y, sigma_z
164 */
165 Scalar U() const { return fU; }
166 Scalar I() const { return fI; }
167 Scalar J() const { return fJ; }
168 Scalar K() const { return fK; }
169
170 // =========== operations ==============
171
172 /**
173 Rotation operation on a cartesian vector
174 */
177
178 const Scalar alpha = fU*fU - fI*fI - fJ*fJ - fK*fK;
179 const Scalar twoQv = 2*(fI*v.X() + fJ*v.Y() + fK*v.Z());
180 const Scalar twoU = 2 * fU;
181 return XYZVector ( alpha * v.X() + twoU * (fJ*v.Z() - fK*v.Y()) + twoQv * fI ,
182 alpha * v.Y() + twoU * (fK*v.X() - fI*v.Z()) + twoQv * fJ ,
183 alpha * v.Z() + twoU * (fI*v.Y() - fJ*v.X()) + twoQv * fK );
184 }
185
186 /**
187 Rotation operation on a displacement vector in any coordinate system
188 */
189 template <class CoordSystem,class Tag>
192 DisplacementVector3D< Cartesian3D<double> > xyz(v.X(), v.Y(), v.Z());
195 vNew.SetXYZ( rxyz.X(), rxyz.Y(), rxyz.Z() );
196 return vNew;
197 }
198
199 /**
200 Rotation operation on a position vector in any coordinate system
201 */
202 template <class CoordSystem, class Tag>
207 return PositionVector3D<CoordSystem,Tag> ( rxyz );
208 }
209
210 /**
211 Rotation operation on a Lorentz vector in any 4D coordinate system
212 */
213 template <class CoordSystem>
217 xyz = operator()(xyz);
218 LorentzVector< PxPyPzE4D<double> > xyzt (xyz.X(), xyz.Y(), xyz.Z(), v.E());
219 return LorentzVector<CoordSystem> ( xyzt );
220 }
221
222 /**
223 Rotation operation on an arbitrary vector v.
224 Preconditions: v must implement methods x(), y(), and z()
225 and the arbitrary vector type must have a constructor taking (x,y,z)
226 */
227 template <class ForeignVector>
228 ForeignVector
229 operator() (const ForeignVector & v) const {
232 return ForeignVector ( rxyz.X(), rxyz.Y(), rxyz.Z() );
233 }
234
235 /**
236 Overload operator * for rotation on a vector
237 */
238 template <class AVector>
239 inline
240 AVector operator* (const AVector & v) const
241 {
242 return operator()(v);
243 }
244
245 /**
246 Invert a rotation in place
247 */
248 void Invert() { fI = -fI; fJ = -fJ; fK = -fK; }
249
250 /**
251 Return inverse of a rotation
252 */
253 Quaternion Inverse() const { return Quaternion(fU, -fI, -fJ, -fK); }
254
255 // ========= Multi-Rotation Operations ===============
256
257 /**
258 Multiply (combine) two rotations
259 */
260 /**
261 Multiply (combine) two rotations
262 */
264 return Quaternion ( fU*q.fU - fI*q.fI - fJ*q.fJ - fK*q.fK ,
265 fU*q.fI + fI*q.fU + fJ*q.fK - fK*q.fJ ,
266 fU*q.fJ - fI*q.fK + fJ*q.fU + fK*q.fI ,
267 fU*q.fK + fI*q.fJ - fJ*q.fI + fK*q.fU );
268 }
269
270 Quaternion operator * (const Rotation3D & r) const;
271 Quaternion operator * (const AxisAngle & a) const;
272 Quaternion operator * (const EulerAngles & e) const;
273 Quaternion operator * (const RotationZYX & r) const;
274 Quaternion operator * (const RotationX & rx) const;
275 Quaternion operator * (const RotationY & ry) const;
276 Quaternion operator * (const RotationZ & rz) const;
277
278 /**
279 Post-Multiply (on right) by another rotation : T = T*R
280 */
281 template <class R>
282 Quaternion & operator *= (const R & r) { return *this = (*this)*r; }
283
284
285 /**
286 Distance between two rotations in Quaternion form
287 Note: The rotation group is isomorphic to a 3-sphere
288 with diametrically opposite points identified.
289 The (rotation group-invariant) is the smaller
290 of the two possible angles between the images of
291 the two totations on that sphere. Thus the distance
292 is never greater than pi/2.
293 */
294
295 Scalar Distance(const Quaternion & q) const ;
296
297 /**
298 Equality/inequality operators
299 */
300 bool operator == (const Quaternion & rhs) const {
301 if( fU != rhs.fU ) return false;
302 if( fI != rhs.fI ) return false;
303 if( fJ != rhs.fJ ) return false;
304 if( fK != rhs.fK ) return false;
305 return true;
306 }
307 bool operator != (const Quaternion & rhs) const {
308 return ! operator==(rhs);
309 }
310
311private:
312
317
318}; // Quaternion
319
320// ============ Class Quaternion ends here ============
321
322/**
323 Distance between two rotations
324 */
325template <class R>
326inline
327typename Quaternion::Scalar
328Distance ( const Quaternion& r1, const R & r2) {return gv_detail::dist(r1,r2);}
329
330/**
331 Multiplication of an axial rotation by an AxisAngle
332 */
333Quaternion operator* (RotationX const & r1, Quaternion const & r2);
334Quaternion operator* (RotationY const & r1, Quaternion const & r2);
335Quaternion operator* (RotationZ const & r1, Quaternion const & r2);
336
337/**
338 Stream Output and Input
339 */
340 // TODO - I/O should be put in the manipulator form
341
342std::ostream & operator<< (std::ostream & os, const Quaternion & q);
343
344
345} // namespace Math
346} // namespace ROOT
347
348#endif // ROOT_Math_GenVector_Quaternion
ROOT::R::TRInterface & r
Definition Object.C:4
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
float * q
typedef void((*Func_t)())
AxisAngle class describing rotation represented with direction axis (3D Vector) and an angle of rotat...
Definition AxisAngle.h:41
DefaultCoordinateSystemTag Default tag for identifying any coordinate system.
Class describing a generic displacement vector in 3 dimensions.
Scalar X() const
Cartesian X, converting if necessary from internal coordinate system.
Scalar Y() const
Cartesian Y, converting if necessary from internal coordinate system.
DisplacementVector3D< CoordSystem, Tag > & SetXYZ(Scalar a, Scalar b, Scalar c)
set the values of the vector from the cartesian components (x,y,z) (if the vector is held in polar or...
Scalar Z() const
Cartesian Z, converting if necessary from internal coordinate system.
EulerAngles class describing rotation as three angles (Euler Angles).
Definition EulerAngles.h:43
Class describing a generic LorentzVector in the 4D space-time, using the specified coordinate system ...
Class describing a generic position vector (point) in 3 dimensions.
Rotation class with the (3D) rotation represented by a unit quaternion (u, i, j, k).
Definition Quaternion.h:47
Scalar Distance(const Quaternion &q) const
Distance between two rotations in Quaternion form Note: The rotation group is isomorphic to a 3-spher...
Quaternion(const OtherRotation &r)
Construct from another supported rotation type (see gv_detail::convert )
Definition Quaternion.h:78
bool operator!=(const Quaternion &rhs) const
Definition Quaternion.h:307
void SetComponents(Scalar u, Scalar i, Scalar j, Scalar k)
Set the components based on four Scalars.
Definition Quaternion.h:149
Quaternion()
Default constructor (identity rotation)
Definition Quaternion.h:58
Quaternion(Scalar u, Scalar i, Scalar j, Scalar k)
Construct from four Scalars representing the coefficients of u, i, j, k.
Definition Quaternion.h:84
DisplacementVector3D< Cartesian3D< double >, DefaultCoordinateSystemTag > XYZVector
Rotation operation on a cartesian vector.
Definition Quaternion.h:175
void Rectify()
Re-adjust components to eliminate small deviations from |Q| = 1 orthonormality.
Quaternion(IT begin, IT end)
Construct given a pair of pointers or iterators defining the beginning and end of an array of four Sc...
Definition Quaternion.h:70
Scalar U() const
Access to the four quaternion components: U() is the coefficient of the identity Pauli matrix,...
Definition Quaternion.h:165
void GetComponents(IT begin, IT end) const
Get the components into data specified by an iterator begin and another to the end of the desired dat...
Definition Quaternion.h:125
Quaternion & operator*=(const R &r)
Post-Multiply (on right) by another rotation : T = T*R.
Definition Quaternion.h:282
Quaternion & operator=(OtherRotation const &r)
Assign from another supported rotation type (see gv_detail::convert )
Definition Quaternion.h:99
bool operator==(const Quaternion &rhs) const
Equality/inequality operators.
Definition Quaternion.h:300
XYZVector operator()(const XYZVector &v) const
Definition Quaternion.h:176
AVector operator*(const AVector &v) const
Overload operator * for rotation on a vector.
Definition Quaternion.h:240
void GetComponents(IT begin) const
Get the components into data specified by an iterator begin.
Definition Quaternion.h:138
Quaternion Inverse() const
Return inverse of a rotation.
Definition Quaternion.h:253
void GetComponents(Scalar &u, Scalar &i, Scalar &j, Scalar &k) const
Get the components into four Scalars.
Definition Quaternion.h:156
void SetComponents(IT begin, IT end)
Set the four components given an iterator to the start of the desired data, and another to the end (4...
Definition Quaternion.h:111
void Invert()
Invert a rotation in place.
Definition Quaternion.h:248
Rotation class with the (3D) rotation represented by a 3x3 orthogonal matrix.
Definition Rotation3D.h:65
Rotation class representing a 3D rotation about the X axis by the angle of rotation.
Definition RotationX.h:43
Rotation class representing a 3D rotation about the Y axis by the angle of rotation.
Definition RotationY.h:43
Rotation class with the (3D) rotation represented by angles describing first a rotation of an angle p...
Definition RotationZYX.h:61
Rotation class representing a 3D rotation about the Z axis by the angle of rotation.
Definition RotationZ.h:43
Namespace for new Math classes and functions.
double dist(Rotation3D const &r1, Rotation3D const &r2)
void convert(R1 const &, R2 const)
AxisAngle operator*(RotationX const &r1, AxisAngle const &r2)
Multiplication of an axial rotation by an AxisAngle.
AxisAngle::Scalar Distance(const AxisAngle &r1, const R &r2)
Distance between two rotations.
Definition AxisAngle.h:320
Rotation3D::Scalar Scalar
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...