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 | Copyright (C) 2011 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
29 {}
30 
32 :
33  t_(t),
34  r_(r)
35 {}
36 
38 :
39  t_(t),
40  r_(quaternion::I)
41 {}
42 
44 :
45  t_(vector::zero),
46  r_(r)
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51 
52 inline const Foam::vector& Foam::septernion::t() const
53 {
54  return t_;
55 }
56 
57 
59 {
60  return r_;
61 }
62 
63 
65 {
66  return t_;
67 }
68 
69 
71 {
72  return r_;
73 }
74 
75 
77 {
78  return t() + r().transform(v);
79 }
80 
81 
83 {
84  return r().invTransform(v - t());
85 }
86 
87 
88 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
89 
91 {
92  t_ = tr.t_;
93  r_ = tr.r_;
94 }
95 
97 {
98  t_ += r().transform(tr.t());
99  r_ *= tr.r();
100 }
101 
102 
103 inline void Foam::septernion::operator=(const vector& t)
104 {
105  t_ = t;
106 }
107 
109 {
110  t_ += t;
111 }
112 
114 {
115  t_ -= t;
116 }
117 
118 
120 {
121  r_ = r;
122 }
123 
125 {
126  r_ *= r;
127 }
128 
130 {
131  r_ /= r;
132 }
133 
134 
135 inline void Foam::septernion::operator*=(const scalar s)
136 {
137  t_ *= s;
138  r_ *= s;
139 }
140 
141 inline void Foam::septernion::operator/=(const scalar s)
142 {
143  t_ /= s;
144  r_ /= s;
145 }
146 
147 
148 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
149 
151 {
152  return septernion(-tr.r().invTransform(tr.t()), conjugate(tr.r()));
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
157 
158 inline bool Foam::operator==(const septernion& tr1, const septernion& tr2)
159 {
160  return (tr1.t() == tr2.t() && tr1.r() == tr2.r());
161 }
162 
163 
164 inline bool Foam::operator!=(const septernion& tr1, const septernion& tr2)
165 {
166  return !operator==(tr1, tr2);
167 }
168 
169 
170 inline Foam::septernion Foam::operator+
171 (
172  const septernion& tr,
173  const vector& t
174 )
175 {
176  return septernion(tr.t() + t, tr.r());
177 }
178 
179 
180 inline Foam::septernion Foam::operator+
181 (
182  const vector& t,
183  const septernion& tr
184 )
185 {
186  return septernion(t + tr.t(), tr.r());
187 }
188 
189 
190 inline Foam::septernion Foam::operator-
191 (
192  const septernion& tr,
193  const vector& t
194 )
195 {
196  return septernion(tr.t() - t, tr.r());
197 }
198 
199 
200 inline Foam::septernion Foam::operator*
201 (
202  const quaternion& r,
203  const septernion& tr
204 )
205 {
206  return septernion(tr.t(), r*tr.r());
207 }
208 
209 
210 inline Foam::septernion Foam::operator*
211 (
212  const septernion& tr,
213  const quaternion& r
214 )
215 {
216  return septernion(tr.t(), tr.r()*r);
217 }
218 
219 
220 inline Foam::septernion Foam::operator/
221 (
222  const septernion& tr,
223  const quaternion& r
224 )
225 {
226  return septernion(tr.t(), tr.r()/r);
227 }
228 
229 
230 inline Foam::septernion Foam::operator*
231 (
232  const septernion& tr1,
233  const septernion& tr2
234 )
235 {
236  return septernion
237  (
238  tr1.t() + tr1.r().transform(tr2.t()),
239  tr1.r().transform(tr2.r())
240  );
241 }
242 
243 
244 inline Foam::septernion Foam::operator/
245 (
246  const septernion& tr1,
247  const septernion& tr2
248 )
249 {
250  return tr1*inv(tr2);
251 }
252 
253 
254 inline Foam::septernion Foam::operator*(const scalar s, const septernion& tr)
255 {
256  return septernion(s*tr.t(), s*tr.r());
257 }
258 
259 
260 inline Foam::septernion Foam::operator*(const septernion& tr, const scalar s)
261 {
262  return septernion(s*tr.t(), s*tr.r());
263 }
264 
265 
266 inline Foam::septernion Foam::operator/(const septernion& tr, const scalar s)
267 {
268  return septernion(tr.t()/s, tr.r()/s);
269 }
270 
271 
272 // ************************************************************************* //
Foam::septernion::operator/=
void operator/=(const quaternion &)
Definition: septernionI.H:129
Foam::septernion
Septernion class used to perform translations and rotations in 3D space.
Definition: septernion.H:64
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::quaternion
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:60
Foam::septernion::operator*=
void operator*=(const septernion &)
Definition: septernionI.H:96
Foam::septernion::transform
vector transform(const vector &v) const
Transform the given vector.
Definition: septernionI.H:76
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:71
Foam::septernion::r
const quaternion & r() const
Definition: septernionI.H:58
Foam::septernion::operator-=
void operator-=(const vector &)
Definition: septernionI.H:113
Foam::I
static const sphericalTensor I(1)
Foam::operator!=
bool operator!=(const particle &, const particle &)
Definition: particle.C:147
Foam::septernion::operator=
void operator=(const septernion &)
Definition: septernionI.H:90
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
Foam::operator*
tmp< fvMatrix< Type > > operator*(const DimensionedField< scalar, volMesh > &, const fvMatrix< Type > &)
Foam::operator/
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
Definition: dimensionedScalar.C:65
Foam::septernion::septernion
septernion()
Construct null.
Definition: septernionI.H:28
Foam::septernion::operator+=
void operator+=(const vector &)
Definition: septernionI.H:108
Foam::septernion::t
const vector & t() const
Definition: septernionI.H:52
Foam::Vector< scalar >
Foam::quaternion::invTransform
vector invTransform(const vector &v) const
Rotate the given vector anti-clockwise.
Definition: quaternionI.H:214
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::conjugate
quaternion conjugate(const quaternion &q)
Return the conjugate of the given quaternion.
Definition: quaternionI.H:359
Foam::septernion::invTransform
vector invTransform(const vector &v) const
Inverse Transform the given vector.
Definition: septernionI.H:82
Foam::quaternion::transform
vector transform(const vector &v) const
Rotate the given vector.
Definition: quaternionI.H:208
Foam::zero
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47