triFaceI.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-2013 OpenFOAM Foundation
9  Copyright (C) 2017-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 #include "IOstreams.H"
30 #include "face.H"
31 #include "triPointRef.H"
32 #include <algorithm> // For std::swap
33 
34 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
35 
36 inline int Foam::triFace::compare(const triFace& a, const triFace& b)
37 {
38  if
39  (
40  (a[0] == b[0] && a[1] == b[1] && a[2] == b[2])
41  || (a[0] == b[1] && a[1] == b[2] && a[2] == b[0])
42  || (a[0] == b[2] && a[1] == b[0] && a[2] == b[1])
43  )
44  {
45  // identical
46  return 1;
47  }
48  else if
49  (
50  (a[0] == b[2] && a[1] == b[1] && a[2] == b[0])
51  || (a[0] == b[1] && a[1] == b[0] && a[2] == b[2])
52  || (a[0] == b[0] && a[1] == b[2] && a[2] == b[1])
53  )
54  {
55  // same face, but reversed orientation
56  return -1;
57  }
58 
59  return 0;
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
66 :
67  FixedList<label, 3>(-1)
68 {}
69 
70 
72 (
73  const label a,
74  const label b,
75  const label c
76 )
77 {
78  operator[](0) = a;
79  operator[](1) = b;
80  operator[](2) = c;
81 }
82 
83 
84 inline Foam::triFace::triFace(std::initializer_list<label> list)
85 :
86  FixedList<label, 3>(list)
87 {}
88 
89 
91 :
92  FixedList<label, 3>(list)
93 {}
94 
95 
97 (
98  const labelUList& list,
99  const FixedList<label, 3>& indices
100 )
101 :
102  FixedList<label, 3>(list, indices)
103 {}
104 
105 
107 :
108  FixedList<label, 3>(is)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 inline Foam::label Foam::triFace::collapse()
115 {
116  // Cannot resize FixedList, so mark duplicates with '-1'
117  // (the lower vertex is retained)
118  // catch any '-1' (eg, if called multiple times)
119 
120  label n = 3;
121  if (operator[](0) == operator[](1) || operator[](1) == -1)
122  {
123  operator[](1) = -1;
124  n--;
125  }
126  else if (operator[](1) == operator[](2) || operator[](2) == -1)
127  {
128  operator[](2) = -1;
129  n--;
130  }
131  if (operator[](0) == operator[](2))
132  {
133  operator[](2) = -1;
134  n--;
135  }
136 
137  return n;
138 }
139 
140 
141 inline void Foam::triFace::flip()
142 {
143  std::swap(operator[](1), operator[](2));
144 }
145 
146 
148 {
149  pointField p(3);
150 
151  p[0] = pts[operator[](0)];
152  p[1] = pts[operator[](1)];
153  p[2] = pts[operator[](2)];
154 
155  return p;
156 }
157 
158 
160 {
161  return Foam::face(*this);
162 }
163 
164 
166 {
167  return triPointRef
168  (
169  points[operator[](0)],
170  points[operator[](1)],
171  points[operator[](2)]
172  );
173 }
174 
175 
177 {
178  return (1.0/3.0)*
179  (
180  points[operator[](0)]
181  + points[operator[](1)]
182  + points[operator[](2)]
183  );
184 }
185 
186 
188 {
189  // As per triPointRef(...).areaNormal()
190  return 0.5*
191  (
192  (points[operator[](1)] - points[operator[](0)])
193  ^(points[operator[](2)] - points[operator[](0)])
194  );
195 }
196 
197 
199 {
200  const vector n(areaNormal(points));
201  const scalar s(Foam::mag(n));
202  return s < ROOTVSMALL ? Zero : n/s;
203 }
204 
205 
206 inline Foam::scalar Foam::triFace::mag(const UList<point>& points) const
207 {
208  return ::Foam::mag(areaNormal(points));
209 }
210 
211 
212 inline Foam::label Foam::triFace::nTriangles() const noexcept
213 {
214  return 1;
215 }
216 
217 
219 {
220  // The starting points of the original and reverse face are identical.
221  return triFace(operator[](0), operator[](2), operator[](1));
222 }
223 
224 
225 inline Foam::label Foam::triFace::which(const label pointLabel) const
226 {
227  return FixedList<label, 3>::find(pointLabel);
228 }
229 
230 
231 inline Foam::label Foam::triFace::thisLabel(const label i) const
232 {
233  return operator[](i);
234 }
235 
236 
237 inline Foam::label Foam::triFace::nextLabel(const label i) const
238 {
239  return operator[]((i == 2 ? 0 : i+1));
240 }
241 
242 
243 inline Foam::label Foam::triFace::prevLabel(const label i) const
244 {
245  return operator[]((i ? i-1 : 2));
246 }
247 
248 
249 inline Foam::scalar Foam::triFace::sweptVol
250 (
251  const UList<point>& opts,
252  const UList<point>& npts
253 ) const
254 {
255  return (1.0/6.0)*
256  (
257  (
258  (npts[operator[](0)] - opts[operator[](0)])
259  & (
260  (opts[operator[](1)] - opts[operator[](0)])
261  ^ (opts[operator[](2)] - opts[operator[](0)])
262  )
263  )
264  + (
265  (npts[operator[](1)] - opts[operator[](1)])
266  & (
267  (opts[operator[](2)] - opts[operator[](1)])
268  ^ (npts[operator[](0)] - opts[operator[](1)])
269  )
270  )
271  + (
272  (opts[operator[](2)] - npts[operator[](2)])
273  & (
274  (npts[operator[](1)] - npts[operator[](2)])
275  ^ (npts[operator[](0)] - npts[operator[](2)])
276  )
277  )
278  );
279 }
280 
281 
283 (
284  const UList<point>& points,
285  const point& refPt,
286  scalar density
287 ) const
288 {
289  // a triangle, do a direct calculation
290  return this->tri(points).inertia(refPt, density);
291 }
292 
293 
295 (
296  const point& p,
297  const vector& q,
298  const UList<point>& points,
299  const intersection::algorithm alg,
300  const intersection::direction dir
301 ) const
302 {
303  return this->tri(points).ray(p, q, alg, dir);
304 }
305 
306 
307 
309 (
310  const point& p,
311  const vector& q,
312  const UList<point>& points,
313  const intersection::algorithm alg,
314  const scalar tol
315 ) const
316 {
317  return this->tri(points).intersection(p, q, alg, tol);
318 }
319 
320 
322 (
323  const point& p,
324  const vector& q,
325  const point& ctr,
326  const UList<point>& points,
327  const intersection::algorithm alg,
328  const scalar tol
329 ) const
330 {
331  return intersection(p, q, points, alg, tol);
332 }
333 
334 
336 (
337  const point& p,
338  const UList<point>& points
339 ) const
340 {
341  return this->tri(points).nearestPoint(p);
342 }
343 
344 
346 (
347  const point& p,
348  const UList<point>& points,
349  label& nearType,
350  label& nearLabel
351 ) const
352 {
353  return this->tri(points).nearestPointClassify(p, nearType, nearLabel);
354 }
355 
356 
357 inline int Foam::triFace::sign
358 (
359  const point& p,
360  const UList<point>& points,
361  const scalar tol
362 ) const
363 {
364  return this->tri(points).sign(p, tol);
365 }
366 
367 
368 inline Foam::label Foam::triFace::nEdges() const noexcept
369 {
370  return 3;
371 }
372 
373 
374 inline Foam::edge Foam::triFace::edge(const label edgei) const
375 {
376  return Foam::edge(thisLabel(edgei), nextLabel(edgei));
377 }
378 
379 
381 (
382  const label edgei,
383  const UList<point>& pts
384 ) const
385 {
386  return vector(pts[nextLabel(edgei)] - pts[thisLabel(edgei)]);
387 }
388 
389 
390 inline Foam::edge Foam::triFace::rcEdge(const label edgei) const
391 {
392  // Edge 0 (forward and reverse) always starts at [0]
393  // for consistency with face flipping
394  const label pointi = edgei ? (3 - edgei) : 0;
395  return Foam::edge(thisLabel(pointi), prevLabel(pointi));
396 }
397 
398 
400 (
401  const label edgei,
402  const UList<point>& pts
403 ) const
404 {
405  // Edge 0 (forward and reverse) always starts at [0]
406  // for consistency with face flipping
407  const label pointi = edgei ? (3 - edgei) : 0;
408  return vector(pts[prevLabel(pointi)] - pts[thisLabel(pointi)]);
409 }
410 
411 
413 {
414  edgeList theEdges(3);
415 
416  theEdges[0].first() = operator[](0);
417  theEdges[0].second() = operator[](1);
418 
419  theEdges[1].first() = operator[](1);
420  theEdges[1].second() = operator[](2);
421 
422  theEdges[2].first() = operator[](2);
423  theEdges[2].second() = operator[](0);
424 
425  return theEdges;
426 }
427 
428 
430 {
431  edgeList theEdges(3);
432 
433  theEdges[0].first() = operator[](0);
434  theEdges[0].second() = operator[](2);
435 
436  theEdges[1].first() = operator[](2);
437  theEdges[1].second() = operator[](1);
438 
439  theEdges[2].first() = operator[](1);
440  theEdges[2].second() = operator[](0);
441 
442  return theEdges;
443 }
444 
445 
446 inline int Foam::triFace::edgeDirection(const Foam::edge& e) const
447 {
448  if (e.first() == operator[](0))
449  {
450  if (e.second() == operator[](1)) return 1; // Forward
451  if (e.second() == operator[](2)) return -1; // Reverse
452  }
453  if (e.first() == operator[](1))
454  {
455  if (e.second() == operator[](2)) return 1; // Forward
456  if (e.second() == operator[](0)) return -1; // Reverse
457  }
458  if (e.first() == operator[](2))
459  {
460  if (e.second() == operator[](0)) return 1; // Forward
461  if (e.second() == operator[](1)) return -1; // Reverse
462  }
463 
464  return 0; // Not found
465 }
466 
467 
468 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
469 
470 inline bool Foam::operator==(const triFace& a, const triFace& b)
471 {
472  return triFace::compare(a,b) != 0;
473 }
474 
475 
476 inline bool Foam::operator!=(const triFace& a, const triFace& b)
477 {
478  return triFace::compare(a,b) == 0;
479 }
480 
481 
482 // ************************************************************************* //
Foam::triFace::inertia
tensor inertia(const UList< point > &points, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triFaceI.H:283
Foam::intersection::direction
direction
Definition: intersection.H:66
Foam::Tensor< scalar >
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::triFace::reverseFace
triFace reverseFace() const
Return face with reverse direction.
Definition: triFaceI.H:218
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Foam::triFace::points
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition: triFaceI.H:147
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::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
Foam::triFace::edge
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition: triFaceI.H:374
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::triFace::thisLabel
label thisLabel(const label i) const
Definition: triFaceI.H:231
Foam::triFace::areaNormal
vector areaNormal(const UList< point > &points) const
The area normal - with magnitude equal to area of face.
Definition: triFaceI.H:187
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
face.H
Foam::triFace::triFaceFace
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:159
Foam::FixedList::find
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: FixedList.C:50
Foam::triFace::unitNormal
vector unitNormal(const UList< point > &points) const
The unit normal.
Definition: triFaceI.H:198
triPointRef.H
Foam::triFace::tri
triPointRef tri(const UList< point > &points) const
Return the triangle.
Definition: triFaceI.H:165
Foam::triFace::flip
void flip()
Flip the face in-place.
Definition: triFaceI.H:141
Foam::triFace::collapse
label collapse()
'Collapse' face by marking duplicate point labels.
Definition: triFaceI.H:114
Foam::triFace::rcEdges
edgeList rcEdges() const
Return list of edges in reverse walk order.
Definition: triFaceI.H:429
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triFace::which
label which(const label pointLabel) const
Find local index on face for the point label, same as find()
Definition: triFaceI.H:225
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::triFace::nTriangles
label nTriangles() const noexcept
Number of triangles after splitting == 1.
Definition: triFaceI.H:212
Foam::triFace::sign
int sign(const point &p, const UList< point > &points, const scalar tol=SMALL) const
The sign for which side of the face plane the point is on.
Definition: triFaceI.H:358
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::triFace::edges
edgeList edges() const
Return list of edges in forward walk order.
Definition: triFaceI.H:412
Foam::intersection::algorithm
algorithm
Definition: intersection.H:72
Foam::triFace::sweptVol
scalar sweptVol(const UList< point > &opts, const UList< point > &npts) const
Return swept-volume from old-points to new-points.
Definition: triFaceI.H:250
Foam::triFace::centre
point centre(const UList< point > &points) const
Return centre (centroid)
Definition: triFaceI.H:176
Foam::triFace::edgeDirection
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition: triFaceI.H:446
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::triFace::nearestPoint
pointHit nearestPoint(const point &p, const UList< point > &points) const
Return nearest point to face.
Definition: triFaceI.H:336
Foam::triFace::nEdges
label nEdges() const noexcept
Return number of edges == 3.
Definition: triFaceI.H:368
Foam::triFace::nearestPointClassify
pointHit nearestPointClassify(const point &p, const UList< point > &points, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Definition: triFaceI.H:346
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
Foam::triFace::rcEdge
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition: triFaceI.H:390
Foam::triFace::ray
pointHit ray(const point &p, const vector &q, const UList< point > &points, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return point intersection with a ray starting at p, in direction q.
Definition: triFaceI.H:295
Foam::Vector< scalar >
Foam::triFace::mag
scalar mag(const UList< point > &points) const
Magnitude of face area.
Definition: triFaceI.H:206
Foam::List< edge >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::triFace::intersection
pointHit intersection(const point &p, const vector &q, const UList< point > &points, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition: triFaceI.H:309
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::intersection
Foam::intersection.
Definition: intersection.H:52
Foam::triFace::prevLabel
label prevLabel(const label i) const
Previous vertex on face.
Definition: triFaceI.H:243
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::triFace::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: triFaceI.H:237
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::triFace::triFace
triFace()
Default construct, with invalid point labels (-1)
Definition: triFaceI.H:65
Foam::triFace::compare
static int compare(const triFace &a, const triFace &b)
Compare triFaces.
Definition: triFaceI.H:36
triFace
face triFace(3)
Foam::triPointRef
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:71