quaternion.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::quaternion
29 
30 Description
31  Quaternion class used to perform rotations in 3D space.
32 
33 SourceFiles
34  quaternionI.H
35  quaternion.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef Foam_quaternion_H
40 #define Foam_quaternion_H
41 
42 #include "scalar.H"
43 #include "vector.H"
44 #include "tensor.H"
45 #include "word.H"
46 #include "Enum.H"
47 #include "stdFoam.H" // For deprecation macros
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 /*---------------------------------------------------------------------------*\
55  Class quaternion Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class quaternion
59 {
60  // Private Data
61 
62  //- Scalar part of the quaternion ( = cos(theta/2) for rotation)
63  scalar w_;
64 
65  //- Vector part of the quaternion ( = axis of rotation)
66  vector v_;
67 
68 
69  //- Multiply vector v by quaternion as if v is a pure quaternion
70  inline quaternion mulq0v(const vector& v) const;
71 
72  //- Conversion of two-axis rotation components into Euler-angles
73  inline static vector twoAxes
74  (
75  const scalar t11,
76  const scalar t12,
77  const scalar c2,
78  const scalar t31,
79  const scalar t32
80  );
81 
82  //- Conversion of three-axis rotation components into Euler-angles
83  inline static vector threeAxes
84  (
85  const scalar t11,
86  const scalar t12,
87  const scalar s2,
88  const scalar t31,
89  const scalar t32
90  );
91 
92 
93 public:
94 
95  // Public Types
96 
97  //- Component type
98  typedef scalar cmptType;
99 
100  //- Magnitude type
101  typedef scalar magType;
102 
103  //- Euler-angle rotation order
104  enum eulerOrder : unsigned char
105  {
106  // Proper Euler angles
107  XZX, XYX, YXY, YZY, ZYZ, ZXZ,
108 
109  // Tait-Bryan angles
111 
112  // Aliases
115  };
116 
117  //- The names for Euler-angle and Tait-Bryan angles,
118  //- including "rollPitchYaw" and "yawPitchRoll" aliases
119  static const Enum<eulerOrder> eulerOrderNames;
120 
122  // Member Constants
123 
124  //- Rank of quaternion is 1
125  static constexpr direction rank = 1;
126 
127 
128  // Static Data Members
129 
130  static constexpr const char* const typeName = "quaternion";
131 
132  static const quaternion zero;
133  static const quaternion I;
134 
135 
136  // Generated Methods
137 
138  //- Default construct
139  quaternion() = default;
141  //- Copy construct
142  quaternion(const quaternion&) = default;
143 
144  //- Copy assignment
145  quaternion& operator=(const quaternion&) = default;
146 
148  // Constructors
149 
150  //- Construct initialized to zero
151  inline quaternion(const Foam::zero);
152 
153  //- Construct given scalar and vector parts
154  inline quaternion(const scalar w, const vector& v);
155 
156  //- Construct rotation quaternion given direction d and angle theta
157  inline quaternion(const vector& d, const scalar theta);
158 
159  //- Construct a rotation quaternion given direction d
160  //- and cosine angle cosTheta and flag if d is normalised
161  inline quaternion
162  (
163  const vector& d,
164  const scalar cosTheta,
165  const bool isNormalised
166  );
167 
168  //- Construct a real quaternion from the given scalar part,
169  //- the vector part = zero
170  inline explicit quaternion(const scalar w);
171 
172  //- Construct a pure imaginary quaternion given the vector part,
173  //- the scalar part = 0
174  inline explicit quaternion(const vector& v);
175 
176  //- Return the unit quaternion (versor) from the given vector
177  //- (w = sqrt(1 - |sqr(v)|))
178  static inline quaternion unit(const vector& v);
179 
180  //- Construct from three Euler rotation angles
181  inline quaternion(const eulerOrder order, const vector& angles);
182 
183  //- Construct from a rotation tensor
184  inline explicit quaternion(const tensor& rotationTensor);
185 
186  //- Construct from Istream
187  explicit quaternion(Istream& is);
188 
189 
190  // Member Functions
191 
192  // Access
193 
194  //- Scalar part of the quaternion ( = cos(theta/2) for rotation)
195  inline scalar w() const noexcept;
196 
197  //- Vector part of the quaternion ( = axis of rotation)
198  inline const vector& v() const noexcept;
199 
200  //- The rotation tensor corresponding to the quaternion
201  inline tensor R() const;
202 
203  //- Return the Euler rotation angles corresponding to the
204  //- specified rotation order
205  inline vector eulerAngles(const eulerOrder order) const;
206 
207 
208  // Edit
209 
210  //- Scalar part of the quaternion ( = cos(theta/2) for rotation)
211  inline scalar& w() noexcept;
212 
213  //- Vector part of the quaternion ( = axis of rotation)
214  inline vector& v() noexcept;
215 
216  //- Inplace normalise the quaternion by its magnitude
217  // For small magnitudes (less than ROOTVSMALL) set to zero.
218  inline quaternion& normalise(const scalar tol = ROOTVSMALL);
219 
220 
221  // Transform
222 
223  //- Rotate the given vector
224  inline vector transform(const vector& v) const;
225 
226  //- Rotate the given vector anti-clockwise
227  inline vector invTransform(const vector& v) const;
228 
229  //- Rotate the given quaternion (and normalise)
230  inline quaternion transform(const quaternion& q) const;
231 
232  //- Rotate the given quaternion anti-clockwise (and normalise)
233  inline quaternion invTransform(const quaternion& q) const;
234 
235 
236  // Member Operators
237 
238  inline void operator+=(const quaternion& q);
239  inline void operator-=(const quaternion& q);
240  inline void operator*=(const quaternion& q);
241  inline void operator/=(const quaternion& q);
242 
243  //- Change scalar portion only
244  inline void operator=(const scalar s);
245 
246  //- Change vector portion only
247  inline void operator=(const vector& v);
248 
249  //- Assign scalar and vector to zero
250  inline void operator=(const Foam::zero);
251 
252  inline void operator*=(const scalar s);
253  inline void operator/=(const scalar s);
254 
255 
256  // Housekeeping
257 
258  //- Inplace normalise the quaternion by its magnitude
259  FOAM_DEPRECATED_FOR(2022-04, "normalise() method")
260  void normalize() { this->normalise(); }
261 
262  //- Return the quaternion normalised by its magnitude
263  FOAM_DEPRECATED_FOR(2022-04, "Use Foam::normalised() function")
264  quaternion normalized() const
265  {
266  quaternion q(*this);
267  q.normalise();
268  return q;
269  }
270 };
271 
272 
273 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
274 
275 //- Contiguous data for quaternion
276 template<> struct is_contiguous<quaternion> : std::true_type {};
277 
278 //- Contiguous scalar data for quaternion
279 template<> struct is_contiguous_scalar<quaternion> : std::true_type {};
280 
281 
282 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
283 
284 inline scalar magSqr(const quaternion& q);
285 inline scalar mag(const quaternion& q);
286 
287 //- Return the conjugate of the given quaternion
288 inline quaternion conjugate(const quaternion& q);
289 
290 //- Return the normalised (unit) quaternion of the given quaternion
291 inline quaternion normalised(const quaternion& q);
292 
293 //- Return the normalised (unit) quaternion of the given quaternion
294 FOAM_DEPRECATED_FOR(2022-04, "Use Foam::normalised() function")
295 inline quaternion normalize(const quaternion& q) { return Foam::normalised(q); }
296 
297 //- Return the inverse of the given quaternion
298 inline quaternion inv(const quaternion& q);
299 
300 //- Return a string representation of a quaternion
301 word name(const quaternion& q);
302 
303 //- Spherical linear interpolation of quaternions
304 quaternion slerp
305 (
306  const quaternion& qa,
307  const quaternion& qb,
308  const scalar t
309 );
310 
311 //- Simple weighted average with sign change
312 quaternion average
313 (
314  const UList<quaternion>& qs,
315  const UList<scalar> w
316 );
317 
318 //- Exponent of a quaternion
319 quaternion exp(const quaternion& q);
320 
321 //- Power of a quaternion
322 quaternion pow(const quaternion& q, const label power);
323 
324 //- Power of a quaternion
325 quaternion pow(const quaternion& q, const scalar power);
326 
327 
328 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
329 
330 Istream& operator>>(Istream& is, quaternion& q);
331 Ostream& operator<<(Ostream& os, const quaternion& q);
333 inline bool operator==(const quaternion& q1, const quaternion& q2);
334 inline bool operator!=(const quaternion& q1, const quaternion& q2);
335 inline quaternion operator+(const quaternion& q1, const quaternion& q2);
336 inline quaternion operator-(const quaternion& q);
337 inline quaternion operator-(const quaternion& q1, const quaternion& q2);
338 inline scalar operator&(const quaternion& q1, const quaternion& q2);
339 inline quaternion operator*(const quaternion& q1, const quaternion& q2);
340 inline quaternion operator/(const quaternion& q1, const quaternion& q2);
341 inline quaternion operator*(const scalar s, const quaternion& q);
342 inline quaternion operator*(const quaternion& q, const scalar s);
343 inline quaternion operator/(const quaternion& q, const scalar s);
344 
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
348 } // End namespace Foam
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 #include "quaternionI.H"
353 
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355 
356 #endif
358 // ************************************************************************* //
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1)
tensor R() const
The rotation tensor corresponding to the quaternion.
Definition: quaternionI.H:352
uint8_t direction
Definition: direction.H:46
quaternion normalized() const
Return the quaternion normalised by its magnitude.
Definition: quaternion.H:338
void normalize()
Inplace normalise the quaternion by its magnitude.
Definition: quaternion.H:332
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector eulerAngles(const eulerOrder order) const
Return the Euler rotation angles corresponding to the specified rotation order.
Definition: quaternionI.H:402
const vector & v() const noexcept
Vector part of the quaternion ( = axis of rotation)
Definition: quaternionI.H:284
scalar cmptType
Component type.
Definition: quaternion.H:105
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static constexpr const char *const typeName
Definition: quaternion.H:145
static quaternion unit(const vector &v)
Return the unit quaternion (versor) from the given vector (w = sqrt(1 - |sqr(v)|)) ...
Definition: quaternionI.H:80
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
tmp< faMatrix< Type > > operator+(const faMatrix< Type > &, const faMatrix< Type > &)
static const Enum< eulerOrder > eulerOrderNames
The names for Euler-angle and Tait-Bryan angles, including "rollPitchYaw" and "yawPitchRoll" aliases...
Definition: quaternion.H:132
class FOAM_DEPRECATED_FOR(2017-05, "Foam::Enum") NamedEnum
Definition: NamedEnum.H:65
static const quaternion zero
Definition: quaternion.H:147
scalar w() const noexcept
Scalar part of the quaternion ( = cos(theta/2) for rotation)
Definition: quaternionI.H:278
quaternion conjugate(const quaternion &q)
Return the conjugate of the given quaternion.
Definition: quaternionI.H:653
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
dimensionedScalar exp(const dimensionedScalar &ds)
tmp< faMatrix< Type > > operator*(const areaScalarField::Internal &, const faMatrix< Type > &)
Istream & operator>>(Istream &, directionInfo &)
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:53
tmp< faMatrix< Type > > operator-(const faMatrix< Type > &)
Unary negation.
eulerOrder
Euler-angle rotation order.
Definition: quaternion.H:115
const direction noexcept
Definition: Scalar.H:258
quaternion & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the quaternion by its magnitude.
Definition: quaternionI.H:302
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
quaternion()=default
Default construct.
OBJstream os(runTime.globalPath()/outputName)
static const quaternion I
Definition: quaternion.H:148
quaternion & operator=(const quaternion &)=default
Copy assignment.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
static constexpr direction rank
Rank of quaternion is 1.
Definition: quaternion.H:140
tmp< GeometricField< Type, faPatchField, areaMesh > > operator &(const faMatrix< Type > &, const DimensionedField< Type, areaMesh > &)
quaternion slerp(const quaternion &qa, const quaternion &qb, const scalar t)
Spherical linear interpolation of quaternions.
Definition: quaternion.C:75
quaternion normalize(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternion.H:379
Includes some standard C++ headers, defines global macros and templates used in multiple places by Op...
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:297
vector invTransform(const vector &v) const
Rotate the given vector anti-clockwise.
Definition: quaternionI.H:331
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Tensor of scalars, i.e. Tensor<scalar>.
scalar magType
Magnitude type.
Definition: quaternion.H:110
vector transform(const vector &v) const
Rotate the given vector.
Definition: quaternionI.H:325
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.