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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 {}
32 
34 :
35  t_(t),
36  r_(r)
37 {}
38 
40 :
41  t_(t),
42  r_(quaternion::I)
43 {}
44 
46 :
47  t_(Zero),
48  r_(r)
49 {}
50 
52 :
53  t_(st.r()),
54  r_(st.E())
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59 
60 inline const Foam::vector& Foam::septernion::t() const
61 {
62  return t_;
63 }
64 
65 
67 {
68  return r_;
69 }
70 
71 
73 {
74  return t_;
75 }
76 
77 
79 {
80  return r_;
81 }
82 
83 
85 {
86  return r().transform(v - t());
87 }
88 
89 
91 {
92  return t() + r().invTransform(v);
93 }
94 
95 
96 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
97 
99 {
100  t_ = tr.t_;
101  r_ = tr.r_;
102 }
103 
105 {
106  t_ = tr.t() + tr.r().invTransform(t_);
107  r_ *= tr.r();
108 }
109 
110 
111 inline void Foam::septernion::operator=(const vector& t)
112 {
113  t_ = t;
114  r_ = quaternion::I;
115 }
116 
118 {
119  t_ += t;
120 }
121 
123 {
124  t_ -= t;
125 }
126 
127 
129 {
130  t_ = Zero;
131  r_ = r;
132 }
133 
135 {
136  t_ = r.invTransform(t_);
137  r_ *= r;
138 }
139 
141 {
142  t_ = r.transform(t_);
143  r_ /= r;
144 }
145 
146 
147 inline void Foam::septernion::operator*=(const scalar s)
148 {
149  t_ *= s;
150  r_ *= s;
151 }
152 
153 inline void Foam::septernion::operator/=(const scalar s)
154 {
155  t_ /= s;
156  r_ /= s;
157 }
158 
159 
160 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
161 
163 {
164  return septernion(-tr.r().transform(tr.t()), conjugate(tr.r()));
165 }
166 
167 
168 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
169 
170 inline bool Foam::operator==(const septernion& tr1, const septernion& tr2)
171 {
172  return (tr1.t() == tr2.t() && tr1.r() == tr2.r());
173 }
174 
175 
176 inline bool Foam::operator!=(const septernion& tr1, const septernion& tr2)
177 {
178  return !operator==(tr1, tr2);
179 }
180 
181 
182 inline Foam::septernion Foam::operator+
183 (
184  const septernion& tr,
185  const vector& t
186 )
187 {
188  return septernion(tr.t() + t, tr.r());
189 }
190 
191 
192 inline Foam::septernion Foam::operator+
193 (
194  const vector& t,
195  const septernion& tr
196 )
197 {
198  return septernion(t + tr.t(), tr.r());
199 }
200 
201 
202 inline Foam::septernion Foam::operator-
203 (
204  const septernion& tr,
205  const vector& t
206 )
207 {
208  return septernion(tr.t() - t, tr.r());
209 }
210 
211 
212 inline Foam::septernion Foam::operator*
213 (
214  const quaternion& r,
215  const septernion& tr
216 )
217 {
218  return septernion(tr.t(), r*tr.r());
219 }
220 
221 
222 inline Foam::septernion Foam::operator*
223 (
224  const septernion& tr,
225  const quaternion& r
226 )
227 {
228  return septernion(r.invTransform(tr.t()), tr.r()*r);
229 }
230 
231 
232 inline Foam::septernion Foam::operator/
233 (
234  const septernion& tr,
235  const quaternion& r
236 )
237 {
238  return septernion(r.transform(tr.t()), tr.r()/r);
239 }
240 
241 
242 inline Foam::septernion Foam::operator*
243 (
244  const septernion& tr1,
245  const septernion& tr2
246 )
247 {
248  return septernion
249  (
250  tr2.r().invTransform(tr1.t()) + tr2.t(),
251  tr1.r().transform(tr2.r())
252  );
253 }
254 
255 
256 inline Foam::septernion Foam::operator/
257 (
258  const septernion& tr1,
259  const septernion& tr2
260 )
261 {
262  return tr1*inv(tr2);
263 }
264 
265 
266 inline Foam::septernion Foam::operator*(const scalar s, const septernion& tr)
267 {
268  return septernion(s*tr.t(), s*tr.r());
269 }
270 
271 
272 inline Foam::septernion Foam::operator*(const septernion& tr, const scalar s)
273 {
274  return septernion(s*tr.t(), s*tr.r());
275 }
276 
277 
278 inline Foam::septernion Foam::operator/(const septernion& tr, const scalar s)
279 {
280  return septernion(tr.t()/s, tr.r()/s);
281 }
282 
283 
284 // ************************************************************************* //
Foam::septernion::operator/=
void operator/=(const quaternion &)
Definition: septernionI.H:140
Foam::septernion
Septernion class used to perform translations and rotations in 3D space.
Definition: septernion.H:66
s
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))
Definition: gmvOutputSpray.H:25
Foam::spatialTransform
Compact representation of the Plücker spatial transformation tensor in terms of the rotation tensor E...
Definition: spatialTransform.H:70
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::quaternion::I
static const quaternion I
Definition: quaternion.H:126
Foam::quaternion
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:56
Foam::septernion::operator*=
void operator*=(const septernion &)
Definition: septernionI.H:104
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:235
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::septernion::t
const vector & t() const
Definition: septernionI.H:60
Foam::septernion::operator-=
void operator-=(const vector &)
Definition: septernionI.H:122
Foam::septernion::operator=
void operator=(const septernion &)
Definition: septernionI.H:98
Foam::operator/
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
Definition: dimensionedScalar.C:68
Foam::septernion::septernion
septernion()
Construct null.
Definition: septernionI.H:30
Foam::septernion::operator+=
void operator+=(const vector &)
Definition: septernionI.H:117
Foam::Vector< scalar >
Foam::operator*
tmp< faMatrix< Type > > operator*(const areaScalarField &, const faMatrix< Type > &)
Foam::septernion::r
const quaternion & r() const
Definition: septernionI.H:66
Foam::quaternion::invTransform
vector invTransform(const vector &v) const
Rotate the given vector anti-clockwise.
Definition: quaternionI.H:337
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:51
Foam::conjugate
quaternion conjugate(const quaternion &q)
Return the conjugate of the given quaternion.
Definition: quaternionI.H:658
Foam::septernion::invTransformPoint
vector invTransformPoint(const vector &v) const
Inverse Transform the given coordinate point.
Definition: septernionI.H:90
Foam::septernion::transformPoint
vector transformPoint(const vector &v) const
Transform the given coordinate point.
Definition: septernionI.H:84
Foam::quaternion::transform
vector transform(const vector &v) const
Rotate the given vector.
Definition: quaternionI.H:331
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95