initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / OpenFOAM / primitives / septernion / septernionI.H
blobeab55ec69370dac41bc3448728e9921d90d4a9f0
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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)
34     t_(t),
35     r_(r)
38 inline Foam::septernion::septernion(const vector& t)
40     t_(t),
41     r_(quaternion::I)
44 inline Foam::septernion::septernion(const quaternion& r)
46     t_(vector::zero),
47     r_(r)
51 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
53 inline const Foam::vector& Foam::septernion::t() const
55     return t_;
59 inline const Foam::quaternion& Foam::septernion::r() const
61     return r_;
65 inline Foam::vector& Foam::septernion::t()
67     return t_;
71 inline Foam::quaternion& Foam::septernion::r()
73     return 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)
93     t_ = tr.t_;
94     r_ = tr.r_;
97 inline void Foam::septernion::operator*=(const septernion& tr)
99     t_ += r().transform(tr.t());
100     r_ *= tr.r();
104 inline void Foam::septernion::operator=(const vector& t)
106     t_ = t;
109 inline void Foam::septernion::operator+=(const vector& t)
111     t_ += t;
114 inline void Foam::septernion::operator-=(const vector& t)
116     t_ -= t;
120 inline void Foam::septernion::operator=(const quaternion& r)
122     r_ = r;
125 inline void Foam::septernion::operator*=(const quaternion& r)
127     r_ *= r;
130 inline void Foam::septernion::operator/=(const quaternion& r)
132     r_ /= r;
136 inline void Foam::septernion::operator*=(const scalar s)
138     t_ *= s;
139     r_ *= s;
142 inline void Foam::septernion::operator/=(const scalar s)
144     t_ /= s;
145     r_ /= 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,
174     const vector& t
177     return septernion(tr.t() + t, tr.r());
181 inline Foam::septernion Foam::operator+
183     const vector& t,
184     const septernion& tr
187     return septernion(t + tr.t(), tr.r());
191 inline Foam::septernion Foam::operator-
193     const septernion& tr,
194     const vector& t
197     return septernion(tr.t() - t, tr.r());
201 inline Foam::septernion Foam::operator*
203     const quaternion& r,
204     const septernion& tr
207     return septernion(tr.t(), r*tr.r());
211 inline Foam::septernion Foam::operator*
213     const septernion& tr,
214     const quaternion& r
217     return septernion(tr.t(), tr.r()*r);
221 inline Foam::septernion Foam::operator/
223     const septernion& tr,
224     const quaternion& r
227     return septernion(tr.t(), tr.r()/r);
231 inline Foam::septernion Foam::operator*
233     const septernion& tr1,
234     const septernion& tr2
237     return septernion
238     (
239         tr1.t() + tr1.r().transform(tr2.t()),
240         tr1.r().transform(tr2.r())
241     );
245 inline Foam::septernion Foam::operator/
247     const septernion& tr1,
248     const septernion& tr2
251     return tr1*inv(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 // ************************************************************************* //