Loading...
Searching...
No Matches
triad.C
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) 2012-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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#include "triad.H"
29#include "transform.H"
30#include "quaternion.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34template<>
35const char* const Foam::triad::vsType::typeName = "triad";
36
37template<>
38const char* const Foam::triad::vsType::componentNames[] = {"x", "y", "z"};
39
40template<>
42(
43 triad::uniform(vector::uniform(0))
44);
45
46template<>
48(
49 triad::uniform(vector::uniform(1))
50);
51
52template<>
54(
55 triad::uniform(vector::uniform(VGREAT))
56);
57
58template<>
60(
61 triad::uniform(vector::uniform(-VGREAT))
62);
63
64template<>
66(
67 triad::uniform(vector::uniform(ROOTVGREAT))
68);
69
70template<>
72(
73 triad::uniform(vector::uniform(-ROOTVGREAT))
74);
75
77(
78 vector(1, 0, 0),
79 vector(0, 1, 0),
80 vector(0, 0, 1)
81);
82
86);
87
88
89// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90
92{
93 tensor Rt(q.R().T());
94 x() = Rt.x();
95 y() = Rt.y();
96 z() = Rt.z();
97}
98
99
100// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101
103{
104 //int which = ((set(0) ? 1 : 0) | (set(1) ? 2 : 0) | (set(2) ? 4 : 0));
105
106 // Hack for 2D z-slab cases
107 // if (!set(2))
108 // {
109 // operator[](2) = vector(0, 0, 1);
110 // }
111
112 // If only two of the axes are set, set the third
113 if (set(0) && set(1) && !set(2))
114 {
115 operator[](2) = orthogonal(operator[](0), operator[](1));
116 }
117 else if (set(0) && set(2) && !set(1))
118 {
119 operator[](1) = orthogonal(operator[](0), operator[](2));
120 }
121 else if (set(1) && set(2) && !set(0))
122 {
123 operator[](0) = orthogonal(operator[](1), operator[](2));
124 }
125
126 // If all the axes are set
127 if (set())
128 {
129 for (int i=0; i<2; i++)
130 {
131 scalar o01 = Foam::mag(operator[](0) & operator[](1));
132 scalar o02 = Foam::mag(operator[](0) & operator[](2));
133 scalar o12 = Foam::mag(operator[](1) & operator[](2));
134
135 if (o01 < o02 && o01 < o12)
136 {
137 operator[](2) = orthogonal(operator[](0), operator[](1));
138
139 // if (o02 < o12)
140 // {
141 // operator[](1) = orthogonal(operator[](0), operator[](2));
142 // }
143 // else
144 // {
145 // operator[](0) = orthogonal(operator[](1), operator[](2));
146 // }
147 }
148 else if (o02 < o12)
149 {
150 operator[](1) = orthogonal(operator[](0), operator[](2));
151
152 // if (o01 < o12)
153 // {
154 // operator[](2) = orthogonal(operator[](0), operator[](1));
155 // }
156 // else
157 // {
158 // operator[](0) = orthogonal(operator[](1), operator[](2));
159 // }
160 }
161 else
162 {
163 operator[](0) = orthogonal(operator[](1), operator[](2));
164
165 // if (o02 < o01)
166 // {
167 // operator[](1) = orthogonal(operator[](0), operator[](2));
168 // }
169 // else
170 // {
171 // operator[](2) = orthogonal(operator[](0), operator[](1));
172 // }
173 }
174 }
175 }
176}
177
178
179void Foam::triad::operator+=(const triad& t2)
180{
181 bool preset[3];
182
183 for (direction i=0; i<3; i++)
184 {
185 if (t2.set(i) && !set(i))
186 {
187 operator[](i) = t2.operator[](i);
188 preset[i] = true;
189 }
190 else
191 {
192 preset[i] = false;
193 }
194 }
195
196 if (set() && t2.set())
197 {
198 direction correspondance[3]{0, 0, 0};
199 short signd[3];
200
201 for (direction i=0; i<3; i++)
202 {
203 if (preset[i])
204 {
205 signd[i] = 0;
206 continue;
207 }
208
209 scalar mostAligned = -1;
210 for (direction j=0; j<3; j++)
211 {
212 bool set = false;
213 for (direction k=0; k<i; k++)
214 {
215 if (correspondance[k] == j)
216 {
217 set = true;
218 break;
219 }
220 }
221
222 if (!set)
223 {
224 scalar a = operator[](i) & t2.operator[](j);
225 scalar maga = Foam::mag(a);
226
227 if (maga > mostAligned)
228 {
229 correspondance[i] = j;
230 mostAligned = maga;
231 signd[i] = sign(a);
232 }
233 }
234 }
236 operator[](i) += signd[i]*t2.operator[](correspondance[i]);
237 }
238 }
239}
240
241
242void Foam::triad::align(const vector& v)
243{
244 if (set())
245 {
246 vector mostAligned
247 (
248 Foam::mag(v & operator[](0)),
249 Foam::mag(v & operator[](1)),
250 Foam::mag(v & operator[](2))
251 );
252
253 scalar mav;
254
255 if
256 (
257 mostAligned.x() > mostAligned.y()
258 && mostAligned.x() > mostAligned.z()
259 )
260 {
261 mav = mostAligned.x();
262 mostAligned = operator[](0);
263 }
264 else if (mostAligned.y() > mostAligned.z())
265 {
266 mav = mostAligned.y();
267 mostAligned = operator[](1);
268 }
269 else
270 {
271 mav = mostAligned.z();
272 mostAligned = operator[](2);
273 }
274
275 if (mav < 0.99)
276 {
277 tensor R(rotationTensor(mostAligned, v));
278
279 operator[](0) = transform(R, operator[](0));
280 operator[](1) = transform(R, operator[](1));
281 operator[](2) = transform(R, operator[](2));
282 }
283 }
284}
285
286
287Foam::triad Foam::triad::sortxyz() const
288{
289 if (!this->set())
290 {
291 return *this;
292 }
293
294 triad t;
295
296 if
297 (
298 Foam::mag(operator[](0).x()) > Foam::mag(operator[](1).x())
299 && Foam::mag(operator[](0).x()) > Foam::mag(operator[](2).x())
300 )
301 {
302 t[0] = operator[](0);
303
304 if (Foam::mag(operator[](1).y()) > Foam::mag(operator[](2).y()))
305 {
306 t[1] = operator[](1);
307 t[2] = operator[](2);
308 }
309 else
310 {
311 t[1] = operator[](2);
312 t[2] = operator[](1);
313 }
314 }
315 else if
316 (
317 Foam::mag(operator[](1).x()) > Foam::mag(operator[](2).x())
318 )
319 {
320 t[0] = operator[](1);
321
322 if (Foam::mag(operator[](0).y()) > Foam::mag(operator[](2).y()))
323 {
324 t[1] = operator[](0);
325 t[2] = operator[](2);
326 }
327 else
328 {
329 t[1] = operator[](2);
330 t[2] = operator[](0);
331 }
332 }
333 else
334 {
335 t[0] = operator[](2);
336
337 if (Foam::mag(operator[](0).y()) > Foam::mag(operator[](1).y()))
338 {
339 t[1] = operator[](0);
340 t[2] = operator[](1);
341 }
342 else
343 {
344 t[1] = operator[](1);
345 t[2] = operator[](0);
346 }
347 }
348
349 if (t[0].x() < 0) t[0] *= -1;
350 if (t[1].y() < 0) t[1] *= -1;
351 if (t[2].z() < 0) t[2] *= -1;
352
353 return t;
354}
355
356
357
358Foam::triad::operator Foam::quaternion() const
359{
360 tensor R;
361
362 R.xx() = x().x();
363 R.xy() = y().x();
364 R.xz() = z().x();
365
366 R.yx() = x().y();
367 R.yy() = y().y();
368 R.yz() = z().y();
369
370 R.zx() = x().z();
371 R.zy() = y().z();
372 R.zz() = z().z();
374 return quaternion(R);
375}
376
377
378// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
379
380Foam::scalar Foam::diff(const triad& A, const triad& B)
381{
382 triad tmpA = A.sortxyz();
383 triad tmpB = B.sortxyz();
384
385 scalar sumDifference = 0;
386
387 for (direction dir = 0; dir < 3; dir++)
388 {
389 if (!tmpA.set(dir) || !tmpB.set(dir))
390 {
391 continue;
392 }
393
394 scalar cosPhi =
395 (tmpA[dir] & tmpB[dir])
396 /(Foam::mag(tmpA[dir])*Foam::mag(tmpA[dir]) + SMALL);
397
398 cosPhi = clamp(cosPhi, -1, 1);
399
400 sumDifference += mag(cosPhi - 1);
401 }
402
403 return (sumDifference/3);
404}
405
406
407// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
scalar y
label k
#define R(A, B, C, D, E, F, K, M)
Tensor< Cmpt > T() const
Return non-Hermitian transpose.
Definition TensorI.H:419
Vector< Cmpt > z() const
Extract vector for row 2.
Definition TensorI.H:160
Vector< Cmpt > y() const
Extract vector for row 1.
Definition TensorI.H:153
Vector< Cmpt > x() const
Extract vector for row 0.
Definition TensorI.H:146
static const char *const componentNames[]
static const Form zero
static const Form one
static const Form max
const Cmpt & operator[](const direction) const
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
static const Form rootMin
static const Form rootMax
friend Ostream & operator(Ostream &, const VectorSpace< Form, Cmpt, Ncmpts > &)
static const Form min
static const char *const typeName
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
const vector & x() const noexcept
Definition Vector.H:135
const vector & z() const noexcept
Definition Vector.H:145
const vector & y() const noexcept
Definition Vector.H:140
Quaternion class used to perform rotations in 3D space.
Definition quaternion.H:54
tensor R() const
The rotation tensor corresponding to the quaternion.
Representation of a 3D Cartesian coordinate system as a Vector of row vectors.
Definition triad.H:63
void align(const vector &v)
Align this triad with the given vector v.
Definition triad.C:235
static const triad unset
Definition triad.H:107
triad()
Default construct as 'unset'.
Definition triadI.H:24
bool set(const direction d) const
Is the vector in the direction d set.
Definition triadI.H:63
triad sortxyz() const
Sort the axes such that they are closest to the x, y and z axes.
Definition triad.C:280
static const triad I
Definition triad.H:106
void operator+=(const triad &t2)
Add the triad t2 to this triad.
Definition triad.C:172
void orthogonalise()
Inplace orthogonalise so that it is ortho-normal.
Definition triad.C:95
static vector orthogonal(const vector &v1, const vector &v2)
Return the vector orthogonal to the two provided.
Definition triadI.H:97
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
T cosPhi(const Vector< T > &a, const Vector< T > &b, const T &tolerance=SMALL)
Calculate angle between a and b in radians.
dimensionedScalar sign(const dimensionedScalar &ds)
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition transform.H:47
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Tensor< scalar > tensor
Definition symmTensor.H:57
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
uint8_t direction
Definition direction.H:49
Vector< scalar > vector
Definition vector.H:57
3D tensor transformation operations.