1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software; you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 inline Foam::septernion::septernion()
32 inline Foam::septernion::septernion(const vector& t, const quaternion& r)
38 inline Foam::septernion::septernion(const vector& t)
44 inline Foam::septernion::septernion(const quaternion& r)
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
53 inline const Foam::vector& Foam::septernion::t() const
59 inline const Foam::quaternion& Foam::septernion::r() const
65 inline Foam::vector& Foam::septernion::t()
71 inline Foam::quaternion& Foam::septernion::r()
77 inline Foam::vector Foam::septernion::transform(const vector& v) const
79 return t() + r().transform(v);
83 inline Foam::vector Foam::septernion::invTransform(const vector& v) const
85 return r().invTransform(v - t());
89 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
91 inline void Foam::septernion::operator=(const septernion& tr)
97 inline void Foam::septernion::operator*=(const septernion& tr)
99 t_ += r().transform(tr.t());
104 inline void Foam::septernion::operator=(const vector& t)
109 inline void Foam::septernion::operator+=(const vector& t)
114 inline void Foam::septernion::operator-=(const vector& t)
120 inline void Foam::septernion::operator=(const quaternion& r)
125 inline void Foam::septernion::operator*=(const quaternion& r)
130 inline void Foam::septernion::operator/=(const quaternion& r)
136 inline void Foam::septernion::operator*=(const scalar s)
142 inline void Foam::septernion::operator/=(const scalar s)
149 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
151 inline Foam::septernion Foam::inv(const septernion& tr)
153 return septernion(-tr.r().invTransform(tr.t()), conjugate(tr.r()));
157 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
159 inline bool Foam::operator==(const septernion& tr1, const septernion& tr2)
161 return (tr1.t() == tr2.t() && tr1.r() == tr2.r());
165 inline bool Foam::operator!=(const septernion& tr1, const septernion& tr2)
167 return !operator==(tr1, tr2);
171 inline Foam::septernion Foam::operator+
173 const septernion& tr,
177 return septernion(tr.t() + t, tr.r());
181 inline Foam::septernion Foam::operator+
187 return septernion(t + tr.t(), tr.r());
191 inline Foam::septernion Foam::operator-
193 const septernion& tr,
197 return septernion(tr.t() - t, tr.r());
201 inline Foam::septernion Foam::operator*
207 return septernion(tr.t(), r*tr.r());
211 inline Foam::septernion Foam::operator*
213 const septernion& tr,
217 return septernion(tr.t(), tr.r()*r);
221 inline Foam::septernion Foam::operator/
223 const septernion& tr,
227 return septernion(tr.t(), tr.r()/r);
231 inline Foam::septernion Foam::operator*
233 const septernion& tr1,
234 const septernion& tr2
239 tr1.t() + tr1.r().transform(tr2.t()),
240 tr1.r().transform(tr2.r())
245 inline Foam::septernion Foam::operator/
247 const septernion& tr1,
248 const septernion& tr2
255 inline Foam::septernion Foam::operator*(const scalar s, const septernion& tr)
257 return septernion(s*tr.t(), s*tr.r());
261 inline Foam::septernion Foam::operator*(const septernion& tr, const scalar s)
263 return septernion(s*tr.t(), s*tr.r());
267 inline Foam::septernion Foam::operator/(const septernion& tr, const scalar s)
269 return septernion(tr.t()/s, tr.r()/s);
273 // ************************************************************************* //