septernionI.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) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30
32:
33 t_(Zero),
34 r_(Zero)
35{}
36
37
39:
40 t_(t),
41 r_(r)
42{}
43
44
46:
47 t_(t),
48 r_(quaternion::I)
49{}
50
51
53:
54 t_(Zero),
55 r_(r)
56{}
57
58
60:
61 t_(st.r()),
62 r_(st.E())
63{}
64
65
66// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67
68inline const Foam::vector& Foam::septernion::t() const
69{
70 return t_;
71}
72
73
75{
76 return r_;
77}
78
79
81{
82 return t_;
83}
84
85
87{
88 return r_;
89}
90
91
93{
94 return r().transform(v - t());
95}
96
97
99{
100 return t() + r().invTransform(v);
101}
102
103
104// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
105
107{
108 t_ = tr.t() + tr.r().invTransform(t_);
109 r_ *= tr.r();
110}
111
112
114{
115 t_ = t;
116 r_ = quaternion::I;
117}
118
119
121{
122 t_ += t;
123}
124
125
127{
128 t_ -= t;
129}
130
131
133{
134 t_ = Zero;
135 r_ = r;
136}
137
138
140{
141 t_ = r.invTransform(t_);
142 r_ *= r;
143}
144
145
147{
148 t_ = r.transform(t_);
149 r_ /= r;
150}
151
152
153inline void Foam::septernion::operator*=(const scalar s)
154{
155 t_ *= s;
156 r_ *= s;
157}
158
159inline void Foam::septernion::operator/=(const scalar s)
160{
161 t_ /= s;
162 r_ /= s;
163}
164
165
166// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
167
169{
170 return septernion(-tr.r().transform(tr.t()), conjugate(tr.r()));
171}
172
173
174// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
175
176inline bool Foam::operator==(const septernion& tr1, const septernion& tr2)
177{
178 return (tr1.t() == tr2.t() && tr1.r() == tr2.r());
179}
180
181
182inline bool Foam::operator!=(const septernion& tr1, const septernion& tr2)
183{
184 return !operator==(tr1, tr2);
185}
186
187
188inline Foam::septernion Foam::operator+
189(
190 const septernion& tr,
191 const vector& t
192)
193{
194 return septernion(tr.t() + t, tr.r());
195}
196
197
198inline Foam::septernion Foam::operator+
199(
200 const vector& t,
201 const septernion& tr
202)
203{
204 return septernion(t + tr.t(), tr.r());
205}
206
207
208inline Foam::septernion Foam::operator-
209(
210 const septernion& tr,
211 const vector& t
212)
213{
214 return septernion(tr.t() - t, tr.r());
215}
216
217
218inline Foam::septernion Foam::operator*
219(
220 const quaternion& r,
221 const septernion& tr
222)
223{
224 return septernion(tr.t(), r*tr.r());
225}
226
227
228inline Foam::septernion Foam::operator*
229(
230 const septernion& tr,
231 const quaternion& r
232)
233{
234 return septernion(r.invTransform(tr.t()), tr.r()*r);
235}
236
237
238inline Foam::septernion Foam::operator/
239(
240 const septernion& tr,
241 const quaternion& r
242)
243{
244 return septernion(r.transform(tr.t()), tr.r()/r);
245}
246
247
248inline Foam::septernion Foam::operator*
249(
250 const septernion& tr1,
251 const septernion& tr2
252)
253{
254 return septernion
255 (
256 tr2.r().invTransform(tr1.t()) + tr2.t(),
257 tr1.r().transform(tr2.r())
258 );
259}
260
261
262inline Foam::septernion Foam::operator/
263(
264 const septernion& tr1,
265 const septernion& tr2
266)
267{
268 return tr1*inv(tr2);
269}
270
271
272inline Foam::septernion Foam::operator*(const scalar s, const septernion& tr)
273{
274 return septernion(s*tr.t(), s*tr.r());
275}
276
277
278inline Foam::septernion Foam::operator*(const septernion& tr, const scalar s)
279{
280 return septernion(s*tr.t(), s*tr.r());
281}
282
283
284inline Foam::septernion Foam::operator/(const septernion& tr, const scalar s)
285{
286 return septernion(tr.t()/s, tr.r()/s);
287}
288
289
290// ************************************************************************* //
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:58
vector transform(const vector &v) const
Rotate the given vector.
Definition: quaternionI.H:331
static const quaternion I
Definition: quaternion.H:132
vector invTransform(const vector &v) const
Rotate the given vector anti-clockwise.
Definition: quaternionI.H:337
Septernion class used to perform translations and rotations in 3D space.
Definition: septernion.H:67
vector invTransformPoint(const vector &v) const
Inverse Transform the given coordinate point.
Definition: septernionI.H:98
void operator*=(const septernion &)
Definition: septernionI.H:106
void operator+=(const vector &)
Definition: septernionI.H:120
septernion()=default
Default construct.
const quaternion & r() const
Definition: septernionI.H:74
vector transformPoint(const vector &v) const
Transform the given coordinate point.
Definition: septernionI.H:92
septernion & operator=(const septernion &)=default
Copy assignment.
void operator-=(const vector &)
Definition: septernionI.H:126
const vector & t() const
Definition: septernionI.H:68
void operator/=(const quaternion &)
Definition: septernionI.H:146
Compact representation of the Plücker spatial transformation tensor in terms of the rotation tensor E...
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
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))
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:94
tmp< faMatrix< Type > > operator*(const areaScalarField::Internal &, const faMatrix< Type > &)
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
quaternion conjugate(const quaternion &q)
Return the conjugate of the given quaternion.
Definition: quaternionI.H:659