triFace.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-2016 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 Class
28  Foam::triFace
29 
30 Description
31  A triangular face using a FixedList of labels corresponding to mesh
32  vertices.
33 
34 See also
35  Foam::face, Foam::triangle
36 
37 SourceFiles
38  triFaceI.H
39  triFaceTemplates.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef Foam_triFace_H
44 #define Foam_triFace_H
45 
46 #include "FixedList.H"
47 #include "edgeList.H"
48 #include "pointHit.H"
49 #include "intersection.H"
50 #include "pointField.H"
51 #include "triangle.H"
52 #include "ListListOps.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 // Forward Declarations
60 class face;
61 class triFace;
62 
63 inline bool operator==(const triFace& a, const triFace& b);
64 inline bool operator!=(const triFace& a, const triFace& b);
65 
66 /*---------------------------------------------------------------------------*\
67  Class triFace Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class triFace
71 :
72  public FixedList<label, 3>
73 {
74 public:
75 
76  // Generated Methods
77 
78  //- The front() accessor (from FixedList) has no purpose
79  void front() = delete;
80 
81  //- The back() accessor (from FixedList) has no purpose
82  void back() = delete;
83 
84 
85  // Constructors
86 
87  //- Default construct, with invalid vertex labels (-1)
88  inline triFace();
89 
90  //- Construct from three vertex labels
91  inline triFace(const label p0, const label p1, const label p2) noexcept;
92 
93  //- Construct from an initializer list of three vertex labels
94  inline explicit triFace(std::initializer_list<label> list);
95 
96  //- Copy construct from a list of three vertex labels.
97  inline explicit triFace(const labelUList& list);
98 
99  //- Copy construct from a subset of vertex labels
100  inline triFace
101  (
102  const labelUList& list,
103  const FixedList<label, 3>& indices
104  );
105 
106  //- Construct from Istream
107  inline triFace(Istream& is);
108 
109 
110  // Member Functions
111 
112  // Access
113 
114  //- The first vertex
115  label a() const noexcept { return get<0>(); }
116 
117  //- The second vertex
118  label b() const noexcept { return get<1>(); }
119 
120  //- The third vertex
121  label c() const noexcept { return get<2>(); }
122 
123  //- The first vertex
124  label& a() noexcept { return get<0>(); }
125 
126  //- The second vertex
127  label& b() noexcept { return get<1>(); }
129  //- The third vertex
130  label& c() noexcept { return get<2>(); }
131 
132 
133  // Queries
134 
135  //- True if vertices are unique and non-negative.
136  inline bool good() const noexcept;
137 
139  // Other
140 
141  //- 'Collapse' face by marking duplicate vertex labels.
142  // Duplicates vertex labels are marked with '-1'
143  // (the lower vertex is retained).
144  // Return the collapsed size.
145  inline label collapse();
146 
147  //- Flip the face in-place.
148  // The starting vertex of the original and flipped face are identical.
149  inline void flip();
150 
151  //- Return the points corresponding to this face
152  inline pointField points(const UList<point>& pts) const;
154  //- Return triangle as a face
155  inline face triFaceFace() const;
156 
157  //- Return the triangle
158  inline triPointRef tri(const UList<point>& points) const;
159 
160  //- Return centre (centroid)
161  inline point centre(const UList<point>& points) const;
162 
163  //- Calculate average value at centroid of face
164  template<class Type>
165  Type average(const UList<point>& unused, const Field<Type>& fld) const;
166 
167  //- The area normal - with magnitude equal to area of face
168  inline vector areaNormal(const UList<point>& points) const;
169 
170  //- The unit normal
171  inline vector unitNormal(const UList<point>& points) const;
172 
173  //- Legacy name for areaNormal()
174  // \deprecated(2018-06) Deprecated for new use
175  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
176  vector normal(const UList<point>& points) const
177  {
178  return areaNormal(points); // Legacy definition
179  }
180 
181  //- Magnitude of face area
182  inline scalar mag(const UList<point>& points) const;
183 
184  //- The enclosing (bounding) box for the face
185  inline Pair<point> box(const UList<point>& points) const;
186 
187  //- Number of triangles after splitting == 1
188  inline label nTriangles() const noexcept;
189 
190  //- Return face with reverse direction
191  // The starting vertex of the original and reverse face are identical.
192  inline triFace reverseFace() const;
193 
194  //- Find local vertex on face for the vertex label, same as find()
195  // \return position in face (0,1,2) or -1 if not found.
196  inline label which(const label vertex) const;
197 
198  //- Next vertex on face
199  inline label nextLabel(const label i) const;
200 
201  //- Previous vertex on face
202  inline label prevLabel(const label i) const;
203 
204  //- The vertex on face - identical to operator[], but with naming
205  //- similar to nextLabel(), prevLabel()
206  inline label thisLabel(const label i) const;
207 
208  //- Return swept-volume from old-points to new-points
209  inline scalar sweptVol
210  (
211  const UList<point>& opts,
212  const UList<point>& npts
213  ) const;
214 
215  //- Return the inertia tensor, with optional reference
216  // point and density specification
217  inline tensor inertia
218  (
219  const UList<point>& points,
220  const point& refPt = vector::zero,
221  scalar density = 1.0
222  ) const;
224  //- Return point intersection with a ray starting at p, in direction q.
225  inline pointHit ray
226  (
227  const point& p,
228  const vector& q,
229  const UList<point>& points,
230  const intersection::algorithm = intersection::FULL_RAY,
231  const intersection::direction dir = intersection::VECTOR
232  ) const;
233 
234  //- Fast intersection with a ray.
235  inline pointHit intersection
236  (
237  const point& p,
238  const vector& q,
239  const UList<point>& points,
240  const intersection::algorithm alg,
241  const scalar tol = 0.0
242  ) const;
243 
244  inline pointHit intersection
245  (
246  const point& p,
247  const vector& q,
248  const point& ctr,
249  const UList<point>& points,
250  const intersection::algorithm alg,
251  const scalar tol = 0.0
252  ) const;
253 
254  //- Return nearest point to face
255  inline pointHit nearestPoint
256  (
257  const point& p,
258  const UList<point>& points
259  ) const;
260 
261 
262  //- Return nearest point to face and classify it:
263  // + near point (nearType=POINT, nearLabel=0, 1, 2)
264  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
265  // Note: edges are counted from starting vertex so
266  // e.g. edge n is from f[n] to f[0], where the face has n + 1
267  // points
269  (
270  const point& p,
271  const UList<point>& points,
272  label& nearType,
273  label& nearLabel
274  ) const;
275 
276  //- The sign for which side of the face plane the point is on.
277  // Uses the supplied tolerance for rounding around zero.
278  // \return
279  // - 0: on plane
280  // - +1: front-side
281  // - -1: back-side
282  inline int sign
283  (
284  const point& p,
285  const UList<point>& points,
286  const scalar tol = SMALL
287  ) const;
288 
289  //- Return number of edges == 3
290  inline label nEdges() const noexcept;
291 
292  //- Return i-th face edge (forward walk order).
293  // The edge 0 is from [0] to [1]
294  // and edge 1 is from [1] to [2]
295  inline Foam::edge edge(const label edgei) const;
296 
297  //- Return vector of i-th face edge (forward walk order).
298  // The edge 0 is from [0] to [1]
299  // and edge 1 is from [1] to [2]
300  inline vector edge(const label edgei, const UList<point>& pts) const;
301 
302  //- Return i-th face edge in reverse walk order.
303  // The rcEdge 0 is from [0] to [n-1]
304  // and rcEdge 1 is from [n-1] to [n-2]
305  inline Foam::edge rcEdge(const label edgei) const;
306 
307  //- Return vector of i-th face edge in reverse walk order.
308  // The rcEdge 0 is from [0] to [n-1]
309  // and rcEdge 1 is from [n-1] to [n-2]
310  inline vector rcEdge(const label edgei, const UList<point>& pts) const;
311 
312  //- Return list of edges in forward walk order.
313  // The edge 0 is from [0] to [1]
314  // and edge 1 is from [1] to [2]
315  inline edgeList edges() const;
316 
317  //- Return list of edges in reverse walk order.
318  // The rcEdge 0 is from [0] to [n-1]
319  // and rcEdge 1 is from [n-1] to [n-2]
320  inline edgeList rcEdges() const;
321 
322 
323  // Searching
324 
325  //- Regular contains(pointi) tests
326  using FixedList<label, 3>::contains;
327 
328  //- True if face contains(edge)
329  inline bool contains(const Foam::edge& e) const;
330 
331  //- Regular find pointi within face
332  using FixedList<label, 3>::find;
333 
334  //- Find the edge within the face.
335  // \return the edge index or -1 if not found
336  inline label find(const Foam::edge& e) const;
337 
338  //- Test the edge direction on the face
339  // \return
340  // - +1: forward (counter-clockwise) on the face
341  // - -1: reverse (clockwise) on the face
342  // - 0: edge not found on the face
343  inline int edgeDirection(const Foam::edge& e) const;
344 
345  //- Compare triFaces
346  // \return:
347  // - 0: different
348  // - +1: identical
349  // - -1: same face, but different orientation
350  static inline int compare(const triFace& a, const triFace& b);
351 
352 
353  // Member Operators
354 
355  //- Increment (offset) vertices by given amount
356  inline void operator+=(const label vertexOffset);
357 
358 
359  // Hashing
360 
361  //- The (commutative) hash value for triFace
362  inline unsigned hash_code(unsigned seed=0) const
363  {
364  // Fortunately we don't need this very often
365  const uLabel t0((*this)[0]);
366  const uLabel t1((*this)[1]);
367  const uLabel t2((*this)[2]);
368 
369  const uLabel val(t0*t1*t2 + t0+t1+t2);
370 
371  return Foam::Hash<uLabel>()(val, seed);
372  }
373 
374  //- Hashing functor for triFace (commutative)
375  struct hasher
376  {
377  unsigned operator()(const triFace& obj, unsigned seed=0) const
378  {
379  return obj.hash_code(seed);
380  }
381  };
382 
383  //- Deprecated(2021-04) hashing functor. Use hasher()
384  // \deprecated(2021-04) - use hasher() functor
385  template<class Unused=bool>
386  struct Hash : triFace::hasher
387  {
388  FOAM_DEPRECATED_FOR(2021-04, "hasher()") Hash() {}
389  };
390 
391 
392  // Housekeeping
393 
394  //- Same as good()
395  bool valid() const noexcept { return good(); }
396 
397  //- Identical to edge()
398  Foam::edge faceEdge(label edgei) const { return this->edge(edgei); }
399 };
400 
401 
402 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
403 
404 //- Contiguous data for triFace
405 template<> struct is_contiguous<triFace> : std::true_type {};
406 
407 //- Contiguous label data for triFace
408 template<> struct is_contiguous_label<triFace> : std::true_type {};
409 
410 //- Hashing for triFace uses commutative (incremental) hash
411 template<> struct Hash<triFace> : triFace::hasher {};
412 
413 
414 //- Specialization to offset faces, used in ListListOps::combineOffset
415 template<>
416 struct offsetOp<triFace>
417 {
418  triFace operator()(const triFace& x, const label offset) const
419  {
420  triFace f(x);
421  f += offset;
422  return f;
423  }
424 };
425 
426 
427 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
428 
429 inline bool operator==(const triFace& a, const triFace& b);
430 inline bool operator!=(const triFace& a, const triFace& b);
431 
432 
433 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
434 
435 } // End namespace Foam
436 
437 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
438 
439 #include "triFaceI.H"
440 
441 #ifdef NoRepository
442  #include "triFaceTemplates.C"
443 #endif
444 
445 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446 
447 #endif
448 
449 // ************************************************************************* //
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
void front()=delete
The front() accessor (from FixedList) has no purpose.
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
uint8_t direction
Definition: direction.H:46
void back()=delete
The back() accessor (from FixedList) has no purpose.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
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
unsigned operator()(const triFace &obj, unsigned seed=0) const
Definition: triFace.H:501
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
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
unsigned hash_code(unsigned seed=0) const
The (commutative) hash value for triFace.
Definition: triFace.H:484
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
bool valid() const noexcept
Same as good()
Definition: triFace.H:524
scalar sweptVol(const UList< point > &opts, const UList< point > &npts) const
Return swept-volume from old-points to new-points.
Definition: triFaceI.H:242
T operator()(const T &x, const label offset) const
Definition: ListListOps.H:101
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 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
Pair< point > box(const UList< point > &points) const
The enclosing (bounding) box for the face.
Definition: triFaceI.H:198
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
FOAM_DEPRECATED_FOR(2021-04, "hasher()") Hash()
Definition: triFace.H:515
const direction noexcept
Definition: Scalar.H:258
label find(const Foam::edge &e) const
Find the edge within the face.
Definition: triFaceI.H:460
Type average(const UList< point > &unused, const Field< Type > &fld) const
Calculate average value at centroid of face.
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
labelList f(nPoints)
vector areaNormal(const UList< point > &points) const
The area normal - with magnitude equal to area of face.
Definition: triFaceI.H:179
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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
Hash function class. The default definition is for primitives. Non-primitives used to hash entries on...
Definition: Hash.H:47
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
void flip()
Flip the face in-place.
Definition: triFaceI.H:143
label c() const noexcept
The third vertex.
Definition: triFace.H:138
bool contains(const Foam::edge &e) const
True if face contains(edge)
Definition: triFaceI.H:482
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
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
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal(const UList< point > &points) const
Legacy name for areaNormal()
Definition: triFace.H:223
Foam::edge faceEdge(label edgei) const
Identical to edge()
Definition: triFace.H:529
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>.
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
Namespace for OpenFOAM.
const pointField & pts