face.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::face
29 
30 Description
31  A face is a list of labels corresponding to mesh vertices.
32 
33 See also
34  Foam::triFace
35 
36 SourceFiles
37  faceI.H
38  face.C
39  faceIntersection.C
40  faceContactSphere.C
41  faceAreaInContact.C
42  faceTemplates.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef Foam_face_H
47 #define Foam_face_H
48 
49 #include "pointField.H"
50 #include "labelList.H"
51 #include "edgeList.H"
52 #include "vectorField.H"
53 #include "faceListFwd.H"
54 #include "intersection.H"
55 #include "pointHit.H"
56 #include "FixedList.H"
57 #include "ListListOps.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 namespace Foam
62 {
63 
64 // Forward Declarations
65 class face;
66 class triFace;
67 template<class T, int SizeMin> class DynamicList;
68 
69 /*---------------------------------------------------------------------------*\
70  Class face Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 class face
74 :
75  public labelList
76 {
77  // Private Member Functions
78 
79  //- Construct list of edge vectors for face
80  tmp<vectorField> calcEdgeVectors
81  (
82  const UList<point>& points
83  ) const;
84 
85  //- Find index of largest internal angle on face
86  label mostConcaveAngle
87  (
88  const UList<point>& points,
89  const vectorField& edges,
90  scalar& edgeCos
91  ) const;
92 
93  //- Enumeration listing the modes for split()
94  enum splitMode
95  {
96  COUNTTRIANGLE,
97  COUNTQUAD,
98  SPLITTRIANGLE,
99  SPLITQUAD
100  };
101 
102  //- Split face into triangles or triangles and quads.
103  // Stores results quadFaces[quadI], triFaces[triI]
104  // \return number of new faces created
105  label split
106  (
107  const splitMode mode,
108  const UList<point>& points,
109  label& triI,
110  label& quadI,
111  faceList& triFaces,
112  faceList& quadFaces
113  ) const;
114 
115 
116 public:
117 
118  // Public Typedefs
119 
120  //- The proximity classifications
121  enum proxType
122  {
123  NONE = 0,
124  POINT,
125  EDGE
126  };
128 
129  // Static Data Members
131  static const char* const typeName;
132 
133 
134  // Constructors
135 
136  //- Default construct
137  constexpr face() noexcept = default;
138 
139  //- Construct given size, with invalid vertex labels (-1)
140  inline explicit face(const label sz);
141 
142  //- Copy construct from list of vertex labels
143  inline explicit face(const labelUList& list);
144 
145  //- Move construct from list of vertex labels
146  inline explicit face(labelList&& list);
147 
148  //- Copy construct from an initializer list of vertex labels
149  inline explicit face(std::initializer_list<label> list);
150 
151  //- Copy construct from list of vertex labels
152  template<unsigned N>
153  inline explicit face(const FixedList<label, N>& list);
154 
155  //- Copy construct from subset of vertex labels
156  inline face(const labelUList& list, const labelUList& indices);
157 
158  //- Copy construct from subset of vertex labels
159  template<unsigned N>
160  inline face(const labelUList& list, const FixedList<label, N>& indices);
161 
162  //- Copy construct from triFace
163  face(const triFace& f);
164 
165  //- Construct from Istream
166  inline face(Istream& is);
167 
168 
169  // Member Functions
170 
171  //- Collapse face by removing duplicate vertex labels
172  // \return the collapsed size
173  label collapse();
174 
175  //- Flip the face in-place.
176  // The starting vertex of the original and flipped face are identical.
177  void flip();
178 
179  //- Return the points corresponding to this face
180  inline pointField points(const UList<point>& pts) const;
181 
182  //- Centre point of face
183  point centre(const UList<point>& points) const;
184 
185  //- Calculate average value at centroid of face
186  template<class Type>
187  Type average
188  (
189  const UList<point>& meshPoints,
190  const Field<Type>& fld
191  ) const;
192 
193  //- The area normal - with magnitude equal to area of face
194  vector areaNormal(const UList<point>& p) const;
195 
196  //- The unit normal
197  inline vector unitNormal(const UList<point>& p) const;
198 
199  //- Legacy name for areaNormal()
200  // \deprecated(2018-06) Deprecated for new use
201  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
202  vector normal(const UList<point>& p) const
203  {
204  return areaNormal(p); // Legacy definition
205  }
206 
207  //- Magnitude of face area
208  inline scalar mag(const UList<point>& p) const;
209 
210  //- The enclosing (bounding) box for the face
211  inline Pair<point> box(const UList<point>& pts) const;
212 
213  //- Return face with reverse direction
214  // The starting vertex of the original and reverse face are identical.
215  face reverseFace() const;
216 
217 
218  // Navigation through face vertices
219 
220  //- Find local vertex on face for the vertex label, same as find()
221  // \return position in face (0,1,2,...) or -1 if not found.
222  inline label which(const label vertex) const;
223 
224  //- The vertex on face - identical to operator[], but with naming
225  //- similar to nextLabel(), prevLabel()
226  inline label thisLabel(const label i) const;
227 
228  //- Next vertex on face
229  inline label nextLabel(const label i) const;
230 
231  //- Previous vertex on face
232  inline label prevLabel(const label i) const;
233 
234 
235  //- Return the volume swept out by the face when its points move
236  scalar sweptVol
237  (
238  const UList<point>& oldPoints,
239  const UList<point>& newPoints
240  ) const;
241 
242  //- Return the inertia tensor, with optional reference
243  //- point and density specification
245  (
246  const UList<point>& p,
247  const point& refPt = vector::zero,
248  scalar density = 1.0
249  ) const;
250 
251  //- Return potential intersection with face with a ray starting
252  //- at p, direction n (does not need to be normalized)
253  // Does face-centre decomposition and returns triangle intersection
254  // point closest to p. Face-centre is calculated from point average.
255  // For a hit, the distance is signed. Positive number
256  // represents the point in front of triangle
257  // In case of miss the point is the nearest point on the face
258  // and the distance is the distance between the intersection point
259  // and the original point.
260  // The half-ray or full-ray intersection and the contact
261  // sphere adjustment of the projection vector is set by the
262  // intersection parameters
263  pointHit ray
264  (
265  const point& p,
266  const vector& n,
267  const UList<point>& meshPoints,
270  ) const;
271 
272  //- Fast intersection with a ray.
273  // Does face-centre decomposition and returns triangle intersection
274  // point closest to p. See triangle::intersection for details.
276  (
277  const point& p,
278  const vector& q,
279  const point& ctr,
280  const UList<point>& meshPoints,
281  const intersection::algorithm alg,
282  const scalar tol = 0.0
283  ) const;
284 
285  //- Return nearest point to face
287  (
288  const point& p,
289  const UList<point>& meshPoints
290  ) const;
291 
292  //- Return nearest point to face and classify it:
293  // + near point (nearType=POINT, nearLabel=0, 1, 2)
294  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
295  // Note: edges are counted from starting vertex so
296  // e.g. edge n is from f[n] to f[0], where the face has n + 1
297  // points
299  (
300  const point& p,
301  const UList<point>& meshPoints,
302  label& nearType,
303  label& nearLabel
304  ) const;
305 
306  //- The sign for the side of the face plane the point is on,
307  //- using three evenly distributed face points for the estimated normal.
308  // Uses the supplied tolerance for rounding around zero.
309  // \return
310  // - 0: on plane
311  // - +1: front-side
312  // - -1: back-side
313  int sign
314  (
315  const point& p,
316  const UList<point>& points,
317  const scalar tol = SMALL
318  ) const;
319 
320  //- Return contact sphere diameter
321  scalar contactSphereDiameter
322  (
323  const point& p,
324  const vector& n,
325  const UList<point>& meshPoints
326  ) const;
327 
328  //- Return area in contact, given the displacement in vertices
329  scalar areaInContact
330  (
331  const UList<point>& meshPoints,
332  const scalarField& v
333  ) const;
334 
335  //- Return number of edges
336  inline label nEdges() const noexcept;
337 
338  //- Return i-th face edge (forward walk order).
339  // The edge 0 is from [0] to [1]
340  // and edge 1 is from [1] to [2]
341  inline Foam::edge edge(const label edgei) const;
342 
343  //- Return vector of i-th face edge (forward walk order).
344  // The edge 0 is from [0] to [1]
345  // and edge 1 is from [1] to [2]
346  inline vector edge(const label edgei, const UList<point>& pts) const;
347 
348  //- Return i-th face edge in reverse walk order.
349  // The rcEdge 0 is from [0] to [n-1]
350  // and rcEdge 1 is from [n-1] to [n-2]
351  inline Foam::edge rcEdge(const label edgei) const;
352 
353  //- Return vector of i-th face edge in reverse walk order.
354  // The rcEdge 0 is from [0] to [n-1]
355  // and rcEdge 1 is from [n-1] to [n-2]
356  inline vector rcEdge(const label edgei, const UList<point>& pts) const;
357 
358  //- Return list of edges in forward walk order.
359  // The edge 0 is from [0] to [1]
360  // and edge 1 is from [1] to [2]
361  edgeList edges() const;
362 
363  //- Return list of edges in reverse walk order.
364  // The rcEdge 0 is from [0] to [n-1]
365  // and rcEdge 1 is from [n-1] to [n-2]
366  edgeList rcEdges() const;
367 
368 
369  // Searching
370 
371  //- Regular contains(pointi) tests
372  using labelList::contains;
373 
374  //- True if face contains(edge)
375  inline bool contains(const Foam::edge& e) const;
376 
377  //- Regular find pointi within face
378  using labelList::find;
379 
380  //- Find the edge within the face.
381  // \return the edge index or -1 if not found
382  label find(const Foam::edge& e) const;
383 
384  //- Test the edge direction on the face
385  // \return
386  // - 0: edge not found on the face
387  // - +1: forward (counter-clockwise) on the face
388  // - -1: reverse (clockwise) on the face
389  int edgeDirection(const Foam::edge& e) const;
390 
391  //- Find the longest edge on a face.
392  label longestEdge(const UList<point>& pts) const;
393 
394 
395  // Face splitting utilities
396 
397  //- Number of triangles after splitting
398  inline label nTriangles() const;
399 
400  //- Number of triangles after splitting
401  label nTriangles(const UList<point>& unused) const;
402 
403  //- Split into triangles using existing points.
404  // Result in triFaces[triI..triI+nTri].
405  // Splits intelligently to maximize triangle quality.
406  // \Return number of faces created.
407  label triangles
408  (
409  const UList<point>& points,
410  label& triI,
411  faceList& triFaces
412  ) const;
413 
414  //- Split into triangles using existing points.
415  // Append to DynamicList.
416  // \Return number of faces created.
417  template<int SizeMin>
418  label triangles
419  (
420  const UList<point>& points,
421  DynamicList<face, SizeMin>& triFaces
422  ) const;
423 
424  //- Number of triangles and quads after splitting
425  // Returns the sum of both
426  label nTrianglesQuads
427  (
428  const UList<point>& points,
429  label& nTris,
430  label& nQuads
431  ) const;
432 
433  //- Split into triangles and quads.
434  // Results in triFaces (starting at triI) and quadFaces
435  // (starting at quadI).
436  // Returns number of new faces created.
437  label trianglesQuads
438  (
439  const UList<point>& points,
440  label& triI,
441  label& quadI,
442  faceList& triFaces,
443  faceList& quadFaces
444  ) const;
445 
446 
447  //- True if the face has at least one vertex in common with other
448  inline bool connected(const labelUList& other) const;
449 
450  //- True if the face has at least one vertex in common with other
451  template<unsigned N>
452  inline bool connected(const FixedList<label, N>& other) const;
453 
454  //- Compare faces
455  // \return
456  // - 0: different
457  // - +1: identical
458  // - -1: same face, but different orientation
459  static int compare(const face& a, const face& b);
460 
461  //- True if the faces have all the same vertices
462  static bool sameVertices(const face& a, const face& b);
463 
464 
465  // Member Operators
466 
467  //- Increment (offset) vertices by given amount
468  inline void operator+=(const label vertexOffset);
469 
470 
471  // Hashing
472 
473  //- The symmetric hash value for face.
474  // Starts with lowest vertex value and walks in the direction
475  // of the next lowest value.
476  static unsigned symmhash_code(const UList<label>& f, unsigned seed=0);
477 
478  //- The hash value for face.
479  // Currently hashes as sequential contiguous data,
480  // but could/should be commutative
481  inline unsigned hash_code(unsigned seed=0) const
482  {
483  return Foam::Hasher(this->cdata(), this->size_bytes(), seed);
484  }
485 
486  //- The symmetric hash value for face.
487  // Starts with lowest vertex value and walks in the direction
488  // of the next lowest value.
489  inline unsigned symmhash_code(unsigned seed=0) const
490  {
491  return face::symmhash_code(*this, seed);
492  }
493 
494  //- Hashing functor for face
495  struct hasher
496  {
497  unsigned operator()(const face& obj, unsigned seed=0) const
498  {
499  return obj.hash_code(seed);
500  }
501  };
502 
503  //- Symmetric hashing functor for face
504  struct symmHasher
505  {
506  unsigned operator()(const face& obj, unsigned seed=0) const
507  {
508  return face::symmhash_code(obj, seed);
509  }
510  };
511 
512 
513  // Housekeeping
514 
515  //- Identical to edge()
516  Foam::edge faceEdge(label edgei) const { return this->edge(edgei); }
517 };
518 
519 
520 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
521 
522 //- Hash face as contiguous (non-commutative) list data
523 template<> struct Hash<face> : face::hasher {};
524 
525 
526 //- Specialization to offset faces, used in ListListOps::combineOffset
527 template<>
528 struct offsetOp<face>
529 {
530  face operator()(const face& x, const label offset) const
531  {
532  face f(x);
533  f += offset;
534  return f;
535  }
536 };
537 
538 
539 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
540 
541 inline bool operator==(const face& a, const face& b);
542 inline bool operator!=(const face& a, const face& b);
543 
544 
545 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
546 
547 } // End namespace Foam
548 
549 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
550 
551 #include "faceI.H"
552 
553 #ifdef NoRepository
554  #include "faceTemplates.C"
555 #endif
556 
557 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
558 
559 #endif
560 
561 // ************************************************************************* //
Close to point.
Definition: face.H:130
scalar areaInContact(const UList< point > &meshPoints, const scalarField &v) const
Return area in contact, given the displacement in vertices.
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition: faceI.H:105
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
proxType
The proximity classifications.
Definition: face.H:127
edgeList edges() const
Return list of edges in forward walk order.
Definition: face.C:764
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
pointHit intersection(const point &p, const vector &q, const point &ctr, const UList< point > &meshPoints, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Unknown proximity.
Definition: face.H:129
label triangles(const UList< point > &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition: face.C:860
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
constexpr face() noexcept=default
Default construct.
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:183
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition: faceI.H:149
vector areaNormal(const UList< point > &p) const
The area normal - with magnitude equal to area of face.
Definition: face.C:569
vector unitNormal(const UList< point > &p) const
The unit normal.
Definition: faceI.H:97
Close to edge.
Definition: face.H:131
Pair< point > box(const UList< point > &pts) const
The enclosing (bounding) box for the face.
Definition: faceI.H:112
static unsigned symmhash_code(const UList< label > &f, unsigned seed=0)
The symmetric hash value for face.
Definition: face.C:414
unsigned operator()(const face &obj, unsigned seed=0) const
Definition: face.H:652
unsigned operator()(const face &obj, unsigned seed=0) const
Definition: face.H:663
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:189
T operator()(const T &x, const label offset) const
Definition: ListListOps.H:101
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal(const UList< point > &p) const
Legacy name for areaNormal()
Definition: face.H:246
label collapse()
Collapse face by removing duplicate vertex labels.
Definition: face.C:467
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition: faceI.H:80
face reverseFace() const
Return face with reverse direction.
Definition: face.C:627
pointHit ray(const point &p, const vector &n, const UList< point > &meshPoints, const intersection::algorithm alg=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return potential intersection with face with a ray starting at p, direction n (does not need to be no...
int sign(const point &p, const UList< point > &points, const scalar tol=SMALL) const
The sign for the side of the face plane the point is on, using three evenly distributed face points f...
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
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
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition: face.C:831
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:65
face triFace(3)
label longestEdge(const UList< point > &pts) const
Find the longest edge on a face.
Definition: face.C:900
scalar sweptVol(const UList< point > &oldPoints, const UList< point > &newPoints) const
Return the volume swept out by the face when its points move.
Definition: face.C:650
label which(const label vertex) const
Find local vertex on face for the vertex label, same as find()
Definition: faceI.H:171
tensor inertia(const UList< point > &p, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference point and density specification.
Definition: face.C:729
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
label trianglesQuads(const UList< point > &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split into triangles and quads.
Definition: face.C:888
label find(const Foam::edge &e) const
Find the edge within the face.
Definition: face.C:806
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
const direction noexcept
Definition: Scalar.H:258
label nEdges() const noexcept
Return number of edges.
Definition: faceI.H:126
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition: faceI.H:133
static const char *const typeName
Definition: face.H:137
pointHit nearestPointClassify(const point &p, const UList< point > &meshPoints, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
pointHit nearestPoint(const point &p, const UList< point > &meshPoints) const
Return nearest point to face.
labelList f(nPoints)
const Vector< label > N(dict.get< Vector< label >>("N"))
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))
Forwards for various types of face lists.
unsigned Hasher(const void *data, size_t len, unsigned seed=0)
Bob Jenkins&#39;s 96-bit mixer hashing function (lookup3)
Definition: Hasher.C:575
Type average(const UList< point > &meshPoints, const Field< Type > &fld) const
Calculate average value at centroid of face.
Definition: faceTemplates.C:46
scalar contactSphereDiameter(const point &p, const vector &n, const UList< point > &meshPoints) const
Return contact sphere diameter.
Foam::edge faceEdge(label edgei) const
Identical to edge()
Definition: face.H:675
bool connected(const labelUList &other) const
True if the face has at least one vertex in common with other.
Definition: faceI.H:201
label nTrianglesQuads(const UList< point > &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition: face.C:874
void flip()
Flip the face in-place.
Definition: face.C:492
bool contains(const Foam::edge &e) const
True if face contains(edge)
Definition: faceI.H:228
edgeList rcEdges() const
Return list of edges in reverse walk order.
Definition: face.C:785
label thisLabel(const label i) const
The vertex on face - identical to operator[], but with naming similar to nextLabel(), prevLabel()
Definition: faceI.H:177
label n
point centre(const UList< point > &points) const
Centre point of face.
Definition: face.C:506
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:297
const label * cdata() const noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:265
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:195
static int compare(const face &a, const face &b)
Compare faces.
Definition: face.C:273
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
static bool sameVertices(const face &a, const face &b)
True if the faces have all the same vertices.
Definition: face.C:374
Tensor of scalars, i.e. Tensor<scalar>.
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: POSIX.C:773
unsigned hash_code(unsigned seed=0) const
The hash value for face.
Definition: face.H:631
std::streamsize size_bytes() const noexcept
Number of contiguous bytes for the List data.
Definition: UListI.H:293
Namespace for OpenFOAM.
const pointField & pts