55 const scalar cosTheta,
56 const bool isNormalised
59 const scalar cosHalfTheta2 = 0.5*(cosTheta + 1);
60 w_ =
sqrt(cosHalfTheta2);
64 v_ =
sqrt(1 - cosHalfTheta2)*d;
95 const eulerOrder order,
111 operator=(quaternion(
vector(0, 0, 1), angles.
x()));
112 operator*=(quaternion(
vector(0, 1, 0), angles.
y()));
113 operator*=(quaternion(
vector(0, 0, 1), angles.
z()));
199 <<
"Unknown euler rotation order "
340 return (
conjugate(*this).mulq0v(u)*(*
this)).v();
361 const scalar
w2 =
sqr(w());
362 const scalar x2 =
sqr(v().
x());
363 const scalar y2 =
sqr(v().
y());
364 const scalar z2 =
sqr(v().z());
366 const scalar txy = 2*v().x()*v().y();
367 const scalar twz = 2*w()*v().z();
368 const scalar txz = 2*v().x()*v().z();
369 const scalar twy = 2*w()*v().y();
370 const scalar tyz = 2*v().y()*v().z();
371 const scalar twx = 2*w()*v().x();
375 w2 + x2 - y2 - z2, txy - twz, txz + twy,
376 txy + twz,
w2 - x2 + y2 - z2, tyz - twx,
377 txz - twy, tyz + twx,
w2 - x2 - y2 + z2
410 const eulerOrder order
413 const scalar
w2 =
sqr(w());
414 const scalar x2 =
sqr(v().
x());
415 const scalar y2 =
sqr(v().
y());
416 const scalar z2 =
sqr(v().z());
424 2*(v().
x()*v().
y() + w()*v().z()),
426 2*(w()*v().
y() - v().
x()*v().z()),
427 2*(v().
y()*v().z() + w()*v().
x()),
437 2*(v().
y()*v().z() - w()*v().
x()),
438 2*(v().
x()*v().z() + w()*v().
y()),
440 2*(v().
y()*v().z() + w()*v().
x()),
441 2*(w()*v().
y() - v().
x()*v().z())
450 2*(w()*v().z() - v().
x()*v().
y()),
452 2*(v().
y()*v().z() + w()*v().
x()),
453 2*(w()*v().
y() - v().
x()*v().z()),
463 2*(v().
x()*v().z() + w()*v().
y()),
464 2*(w()*v().
x() - v().
y()*v().z()),
466 2*(v().
x()*v().z() - w()*v().
y()),
467 2*(v().
y()*v().z() + w()*v().
x())
476 2*(v().
x()*v().z() + w()*v().
y()),
478 2*(w()*v().
x() - v().
y()*v().z()),
479 2*(v().
x()*v().
y() + w()*v().z()),
489 2*(v().
x()*v().
y() - w()*v().z()),
490 2*(v().
y()*v().z() + w()*v().
x()),
492 2*(v().
x()*v().
y() + w()*v().z()),
493 2*(w()*v().
x() - v().
y()*v().z())
502 2*(w()*v().
y() - v().
x()*v().z()),
504 2*(v().
x()*v().
y() + w()*v().z()),
505 2*(w()*v().
x() - v().
y()*v().z()),
515 2*(v().
y()*v().z() + w()*v().
x()),
516 2*(w()*v().z() - v().
x()*v().
y()),
518 2*(v().
y()*v().z() - w()*v().
x()),
519 2*(v().
x()*v().
y() + w()*v().z())
528 2*(w()*v().
x() - v().
y()*v().z()),
530 2*(v().
x()*v().z() + w()*v().
y()),
531 2*(w()*v().z() - v().
x()*v().
y()),
541 2*(v().
x()*v().
y() + w()*v().z()),
542 2*(w()*v().
y() - v().
x()*v().z()),
544 2*(v().
x()*v().
y() - w()*v().z()),
545 2*(v().
x()*v().z() + w()*v().
y())
554 2*(v().
y()*v().z() + w()*v().
x()),
556 2*(w()*v().z() - v().
x()*v().
y()),
557 2*(v().
x()*v().z() + w()*v().
y()),
567 2*(v().
x()*v().z() - w()*v().
y()),
568 2*(v().
x()*v().
y() + w()*v().z()),
570 2*(v().
x()*v().z() + w()*v().
y()),
571 2*(w()*v().z() - v().
x()*v().
y())
578 <<
"Unknown euler rotation order "
604 w() =
w()*q.
w() - (
v() & q.
v());
605 v() =
w0*q.
v() + q.
w()*
v() + (
v() ^ q.
v());
712 const quaternion& q1,
716 return quaternion(q1.w() + q2.w(), q1.v() + q2.v());
728 const quaternion& q1,
732 return quaternion(q1.w() - q2.w(), q1.v() - q2.v());
738 return q1.
w()*q2.
w() + (q1.
v() & q2.
v());
744 const quaternion& q1,
750 q1.w()*q2.w() - (q1.v() & q2.v()),
751 q1.w()*q2.v() + q2.w()*q1.v() + (q1.v() ^ q2.v())
758 const quaternion& q1,
780 return quaternion(q.w()/
s, q.v()/
s);
const Cmpt & yy() const noexcept
const Cmpt & zy() const noexcept
const Cmpt & yx() const noexcept
const Cmpt & xx() const noexcept
const Cmpt & yz() const noexcept
const Cmpt & zz() const noexcept
const Cmpt & xy() const noexcept
const Cmpt & xz() const noexcept
const Cmpt & zx() const noexcept
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
Quaternion class used to perform rotations in 3D space.
void operator/=(const quaternion &q)
quaternion & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the quaternion by its magnitude.
vector eulerAngles(const eulerOrder order) const
Return the Euler rotation angles corresponding to the specified rotation order.
const vector & v() const noexcept
Vector part of the quaternion ( = axis of rotation).
tensor R() const
The rotation tensor corresponding to the quaternion.
void operator-=(const quaternion &q)
static quaternion unit(const vector &v)
Return the unit quaternion (versor) from the given vector (w = sqrt(1 - |sqr(v)|)).
void operator+=(const quaternion &q)
void operator*=(const quaternion &q)
scalar w() const noexcept
Scalar part of the quaternion ( = cos(theta/2) for rotation).
quaternion & operator=(const quaternion &)=default
Copy assignment.
vector transform(const vector &v) const
Rotate the given vector.
quaternion()=default
Default construct.
eulerOrder
Euler-angle rotation order.
vector invTransform(const vector &v) const
Rotate the given vector anti-clockwise.
Tensor of scalars, i.e. Tensor<scalar>.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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))
tmp< faMatrix< Type > > operator-(const faMatrix< Type > &)
Unary negation.
bool operator!=(const eddy &a, const eddy &b)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
dimensionedScalar sin(const dimensionedScalar &ds)
bool equal(const T &a, const T &b)
Compare two values for equality.
tmp< faMatrix< Type > > operator*(const areaScalarField::Internal &, const faMatrix< Type > &)
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< GeometricField< Type, faPatchField, areaMesh > > operator&(const faMatrix< Type > &, const DimensionedField< Type, areaMesh > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
quaternion conjugate(const quaternion &q)
Return the conjugate of the given quaternion.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)