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
84(
85 triad::uniform(vector::uniform(VGREAT))
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 // Hack for 2D z-slab cases
105 // if (!set(2))
106 // {
107 // operator[](2) = vector(0, 0, 1);
108 // }
109
110 // If only two of the axes are set, set the third
111 if (set(0) && set(1) && !set(2))
112 {
113 operator[](2) = orthogonal(operator[](0), operator[](1));
114 }
115 else if (set(0) && set(2) && !set(1))
116 {
117 operator[](1) = orthogonal(operator[](0), operator[](2));
118 }
119 else if (set(1) && set(2) && !set(0))
120 {
121 operator[](0) = orthogonal(operator[](1), operator[](2));
122 }
123
124 // If all the axes are set
125 if (set())
126 {
127 for (int i=0; i<2; i++)
128 {
129 scalar o01 = mag(operator[](0) & operator[](1));
130 scalar o02 = mag(operator[](0) & operator[](2));
131 scalar o12 = mag(operator[](1) & operator[](2));
132
133 if (o01 < o02 && o01 < o12)
134 {
135 operator[](2) = orthogonal(operator[](0), operator[](1));
136
137 // if (o02 < o12)
138 // {
139 // operator[](1) = orthogonal(operator[](0), operator[](2));
140 // }
141 // else
142 // {
143 // operator[](0) = orthogonal(operator[](1), operator[](2));
144 // }
145 }
146 else if (o02 < o12)
147 {
148 operator[](1) = orthogonal(operator[](0), operator[](2));
149
150 // if (o01 < o12)
151 // {
152 // operator[](2) = orthogonal(operator[](0), operator[](1));
153 // }
154 // else
155 // {
156 // operator[](0) = orthogonal(operator[](1), operator[](2));
157 // }
158 }
159 else
160 {
161 operator[](0) = orthogonal(operator[](1), operator[](2));
162
163 // if (o02 < o01)
164 // {
165 // operator[](1) = orthogonal(operator[](0), operator[](2));
166 // }
167 // else
168 // {
169 // operator[](2) = orthogonal(operator[](0), operator[](1));
170 // }
171 }
172 }
173 }
174}
175
176
178{
179 bool preset[3];
180
181 for (direction i=0; i<3; i++)
182 {
183 if (t2.set(i) && !set(i))
184 {
185 operator[](i) = t2.operator[](i);
186 preset[i] = true;
187 }
188 else
189 {
190 preset[i] = false;
191 }
192 }
193
194 if (set() && t2.set())
195 {
196 direction correspondance[3]{0, 0, 0};
197 short signd[3];
198
199 for (direction i=0; i<3; i++)
200 {
201 if (preset[i])
202 {
203 signd[i] = 0;
204 continue;
205 }
206
207 scalar mostAligned = -1;
208 for (direction j=0; j<3; j++)
209 {
210 bool set = false;
211 for (direction k=0; k<i; k++)
212 {
213 if (correspondance[k] == j)
214 {
215 set = true;
216 break;
217 }
218 }
219
220 if (!set)
221 {
222 scalar a = operator[](i) & t2.operator[](j);
223 scalar maga = mag(a);
224
225 if (maga > mostAligned)
226 {
227 correspondance[i] = j;
228 mostAligned = maga;
229 signd[i] = sign(a);
230 }
231 }
232 }
233
234 operator[](i) += signd[i]*t2.operator[](correspondance[i]);
235 }
236 }
237}
238
239
241{
242 if (set())
243 {
244 vector mostAligned
245 (
246 mag(v & operator[](0)),
247 mag(v & operator[](1)),
248 mag(v & operator[](2))
249 );
250
251 scalar mav;
252
253 if
254 (
255 mostAligned.x() > mostAligned.y()
256 && mostAligned.x() > mostAligned.z()
257 )
258 {
259 mav = mostAligned.x();
260 mostAligned = operator[](0);
261 }
262 else if (mostAligned.y() > mostAligned.z())
263 {
264 mav = mostAligned.y();
265 mostAligned = operator[](1);
266 }
267 else
268 {
269 mav = mostAligned.z();
270 mostAligned = operator[](2);
271 }
272
273 if (mav < 0.99)
274 {
275 tensor R(rotationTensor(mostAligned, v));
276
277 operator[](0) = transform(R, operator[](0));
278 operator[](1) = transform(R, operator[](1));
279 operator[](2) = transform(R, operator[](2));
280 }
281 }
282}
283
284
286{
287 if (!this->set())
288 {
289 return *this;
290 }
291
292 triad t;
293
294 if
295 (
296 mag(operator[](0).x()) > mag(operator[](1).x())
297 && mag(operator[](0).x()) > mag(operator[](2).x())
298 )
299 {
300 t[0] = operator[](0);
301
302 if (mag(operator[](1).y()) > mag(operator[](2).y()))
303 {
304 t[1] = operator[](1);
305 t[2] = operator[](2);
306 }
307 else
308 {
309 t[1] = operator[](2);
310 t[2] = operator[](1);
311 }
312 }
313 else if
314 (
315 mag(operator[](1).x()) > mag(operator[](2).x())
316 )
317 {
318 t[0] = operator[](1);
319
320 if (mag(operator[](0).y()) > mag(operator[](2).y()))
321 {
322 t[1] = operator[](0);
323 t[2] = operator[](2);
324 }
325 else
326 {
327 t[1] = operator[](2);
328 t[2] = operator[](0);
329 }
330 }
331 else
332 {
333 t[0] = operator[](2);
334
335 if (mag(operator[](0).y()) > mag(operator[](1).y()))
336 {
337 t[1] = operator[](0);
338 t[2] = operator[](1);
339 }
340 else
341 {
342 t[1] = operator[](1);
343 t[2] = operator[](0);
344 }
345 }
346
347 if (t[0].x() < 0) t[0] *= -1;
348 if (t[1].y() < 0) t[1] *= -1;
349 if (t[2].z() < 0) t[2] *= -1;
350
351 return t;
352}
353
354
355
357{
358 tensor R;
359
360 R.xx() = x().x();
361 R.xy() = y().x();
362 R.xz() = z().x();
363
364 R.yx() = x().y();
365 R.yy() = y().y();
366 R.yz() = z().y();
367
368 R.zx() = x().z();
369 R.zy() = y().z();
370 R.zz() = z().z();
371
372 return quaternion(R);
373}
374
375
376// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
377
378Foam::scalar Foam::diff(const triad& A, const triad& B)
379{
380 triad tmpA = A.sortxyz();
381 triad tmpB = B.sortxyz();
382
383 scalar sumDifference = 0;
384
385 for (direction dir = 0; dir < 3; dir++)
386 {
387 if (!tmpA.set(dir) || !tmpB.set(dir))
388 {
389 continue;
390 }
391
392 scalar cosPhi =
393 (tmpA[dir] & tmpB[dir])
394 /(mag(tmpA[dir])*mag(tmpA[dir]) + SMALL);
395
396 cosPhi = min(max(cosPhi, -1), 1);
397
398 sumDifference += mag(cosPhi - 1);
399 }
400
401 return (sumDifference/3);
402}
403
404
405// ************************************************************************* //
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)
Vector< Cmpt > z() const
Extract vector for row 2.
Definition: TensorI.H:293
Vector< Cmpt > y() const
Extract vector for row 1.
Definition: TensorI.H:286
Vector< Cmpt > x() const
Extract vector for row 0.
Definition: TensorI.H:279
static const Form zero
Definition: VectorSpace.H:115
static const char *const componentNames[]
Definition: VectorSpace.H:114
static const Form one
Definition: VectorSpace.H:116
static const Form rootMax
Definition: VectorSpace.H:119
static const Form max
Definition: VectorSpace.H:117
static const Form rootMin
Definition: VectorSpace.H:120
static const Form min
Definition: VectorSpace.H:118
static const char *const typeName
Definition: VectorSpace.H:113
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
const vector & z() const
Access to the vector z component.
Definition: VectorI.H:85
const vector & y() const
Access to the vector y component.
Definition: VectorI.H:79
const vector & x() const
Access to the vector x component.
Definition: VectorI.H:73
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:58
tensor R() const
The rotation tensor corresponding to the quaternion.
Definition: quaternionI.H:358
Representation of a 3D Cartesian coordinate system as a Vector of row vectors.
Definition: triad.H:67
void align(const vector &v)
Align this triad with the given vector v.
Definition: triad.C:240
static const triad unset
Definition: triad.H:97
triad()
Default construct.
Definition: triadI.H:31
bool set(const direction d) const
Is the vector in the direction d set.
Definition: triadI.H:70
triad sortxyz() const
Sort the axes such that they are closest to the x, y and z axes.
Definition: triad.C:285
static const triad I
Definition: triad.H:96
void operator+=(const triad &t2)
Add the triad t2 to this triad.
Definition: triad.C:177
void orthogonalise()
Inplace orthogonalise so that it is ortho-normal.
Definition: triad.C:102
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
dimensionedScalar sign(const dimensionedScalar &ds)
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:51
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:378
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
uint8_t direction
Definition: direction.H:56
3D tensor transformation operations.