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