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-2020 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 face_H
47 #define 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  //- Edge to the right of face vertex i
80  inline label right(const label i) const;
81 
82  //- Edge to the left of face vertex i
83  inline label left(const label i) const;
84 
85  //- Construct list of edge vectors for face
86  tmp<vectorField> calcEdges
87  (
88  const UList<point>& points
89  ) const;
90 
91  //- Cos between neighbouring edges
92  scalar edgeCos
93  (
94  const vectorField& edges,
95  const label index
96  ) const;
97 
98  //- Find index of largest internal angle on face
99  label mostConcaveAngle
100  (
101  const UList<point>& points,
102  const vectorField& edges,
103  scalar& edgeCos
104  ) const;
105 
106  //- Enumeration listing the modes for split()
107  enum splitMode
108  {
109  COUNTTRIANGLE,
110  COUNTQUAD,
111  SPLITTRIANGLE,
112  SPLITQUAD
113  };
114 
115  //- Split face into triangles or triangles and quads.
116  // Stores results quadFaces[quadI], triFaces[triI]
117  // \return number of new faces created
118  label split
119  (
120  const splitMode mode,
121  const UList<point>& points,
122  label& triI,
123  label& quadI,
124  faceList& triFaces,
125  faceList& quadFaces
126  ) const;
127 
128 
129 public:
130 
131  // Public Typedefs
132 
133  //- The proximity classifications
134  enum proxType
135  {
136  NONE = 0,
138  EDGE
139  };
140 
141 
142  // Static Data Members
143 
144  static const char* const typeName;
145 
146 
147  // Constructors
148 
149  //- Default construct
150  constexpr face() noexcept = default;
151 
152  //- Construct given size, with invalid point labels (-1)
153  inline explicit face(const label sz);
154 
155  //- Copy construct from list of point labels
156  inline explicit face(const labelUList& list);
157 
158  //- Move construct from list of point labels
159  inline explicit face(labelList&& list);
160 
161  //- Copy construct from an initializer list of point labels
162  inline explicit face(std::initializer_list<label> list);
163 
164  //- Copy construct from list of point labels
165  template<unsigned N>
166  inline explicit face(const FixedList<label, N>& list);
167 
168  //- Copy construct from subset of point labels
169  inline face(const labelUList& list, const labelUList& indices);
170 
171  //- Copy construct from subset of point labels
172  template<unsigned N>
173  inline face(const labelUList& list, const FixedList<label, N>& indices);
174 
175  //- Copy construct from triFace
176  face(const triFace& f);
177 
178  //- Construct from Istream
179  inline face(Istream& is);
180 
181 
182  // Member Functions
183 
184  //- Collapse face by removing duplicate point labels
185  // \return the collapsed size
186  label collapse();
187 
188  //- Flip the face in-place.
189  // The starting points of the original and flipped face are identical.
190  void flip();
191 
192  //- Return the points corresponding to this face
193  inline pointField points(const UList<point>& points) const;
194 
195  //- Centre point of face
196  point centre(const UList<point>& points) const;
197 
198  //- Calculate average value at centroid of face
199  template<class Type>
200  Type average
201  (
202  const UList<point>& meshPoints,
203  const Field<Type>& fld
204  ) const;
205 
206  //- The area normal - with magnitude equal to area of face
207  vector areaNormal(const UList<point>& p) const;
208 
209  //- The unit normal
210  inline vector unitNormal(const UList<point>& p) const;
211 
212  //- Legacy name for areaNormal()
213  // \deprecated(2018-06) Deprecated for new use
214  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
215  vector normal(const UList<point>& p) const
216  {
217  return areaNormal(p); // Legacy definition
218  }
219 
220  //- Magnitude of face area
221  inline scalar mag(const UList<point>& p) const;
222 
223  //- Return face with reverse direction
224  // The starting points of the original and reverse face are identical.
225  face reverseFace() const;
226 
227 
228  // Navigation through face vertices
229 
230  //- Return true if the point label is found in face.
231  inline bool found(const label pointLabel) const;
232 
233  //- Find local index on face for the point label,
234  // \return position in face (0,1,2,...) or -1 if not found.
235  inline label which(const label pointLabel) const;
236 
237  //- Next vertex on face
238  inline label nextLabel(const label i) const;
239 
240  //- Previous vertex on face
241  inline label prevLabel(const label i) const;
242 
243 
244  //- Return the volume swept out by the face when its points move
245  scalar sweptVol
246  (
247  const UList<point>& oldPoints,
248  const UList<point>& newPoints
249  ) const;
250 
251  //- Return the inertia tensor, with optional reference
252  //- point and density specification
254  (
255  const UList<point>& p,
256  const point& refPt = vector::zero,
257  scalar density = 1.0
258  ) const;
259 
260  //- Return potential intersection with face with a ray starting
261  //- at p, direction n (does not need to be normalized)
262  // Does face-centre decomposition and returns triangle intersection
263  // point closest to p. Face-centre is calculated from point average.
264  // For a hit, the distance is signed. Positive number
265  // represents the point in front of triangle
266  // In case of miss the point is the nearest point on the face
267  // and the distance is the distance between the intersection point
268  // and the original point.
269  // The half-ray or full-ray intersection and the contact
270  // sphere adjustment of the projection vector is set by the
271  // intersection parameters
272  pointHit ray
273  (
274  const point& p,
275  const vector& n,
276  const UList<point>& meshPoints,
279  ) const;
280 
281  //- Fast intersection with a ray.
282  // Does face-centre decomposition and returns triangle intersection
283  // point closest to p. See triangle::intersection for details.
285  (
286  const point& p,
287  const vector& q,
288  const point& ctr,
289  const UList<point>& meshPoints,
290  const intersection::algorithm alg,
291  const scalar tol = 0.0
292  ) const;
293 
294  //- Return nearest point to face
296  (
297  const point& p,
298  const UList<point>& meshPoints
299  ) const;
300 
301  //- Return nearest point to face and classify it:
302  // + near point (nearType=POINT, nearLabel=0, 1, 2)
303  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
304  // Note: edges are counted from starting vertex so
305  // e.g. edge n is from f[n] to f[0], where the face has n + 1
306  // points
308  (
309  const point& p,
310  const UList<point>& meshPoints,
311  label& nearType,
312  label& nearLabel
313  ) const;
314 
315  //- The sign for the side of the face plane the point is on,
316  //- using three evenly distributed face points for the estimated normal.
317  // Uses the supplied tolerance for rounding around zero.
318  // \return
319  // - 0: on plane
320  // - +1: front-side
321  // - -1: back-side
322  int sign
323  (
324  const point& p,
325  const UList<point>& points,
326  const scalar tol = SMALL
327  ) const;
328 
329  //- Return contact sphere diameter
330  scalar contactSphereDiameter
331  (
332  const point& p,
333  const vector& n,
334  const UList<point>& meshPoints
335  ) const;
336 
337  //- Return area in contact, given the displacement in vertices
338  scalar areaInContact
339  (
340  const UList<point>& meshPoints,
341  const scalarField& v
342  ) const;
343 
344  //- Return number of edges
345  inline label nEdges() const;
346 
347  //- Return edges in face point ordering,
348  // i.e. edges()[0] is edge between [0] and [1]
349  edgeList edges() const;
350 
351  //- Return n-th face edge
352  inline edge faceEdge(const label n) const;
353 
354  //- The edge direction on the face
355  // \return
356  // - 0: edge not found on the face
357  // - +1: forward (counter-clockwise) on the face
358  // - -1: reverse (clockwise) on the face
359  int edgeDirection(const edge& e) const;
360 
361  //- Find the longest edge on a face.
362  label longestEdge(const UList<point>& pts) const;
363 
364  // Face splitting utilities
365 
366  //- Number of triangles after splitting
367  inline label nTriangles() const;
368 
369  //- Number of triangles after splitting
370  label nTriangles(const UList<point>& unused) const;
371 
372  //- Split into triangles using existing points.
373  // Result in triFaces[triI..triI+nTri].
374  // Splits intelligently to maximize triangle quality.
375  // \Return number of faces created.
376  label triangles
377  (
378  const UList<point>& points,
379  label& triI,
380  faceList& triFaces
381  ) const;
382 
383  //- Split into triangles using existing points.
384  // Append to DynamicList.
385  // \Return number of faces created.
386  template<int SizeMin>
387  label triangles
388  (
389  const UList<point>& points,
391  ) const;
392 
393  //- Number of triangles and quads after splitting
394  // Returns the sum of both
395  label nTrianglesQuads
396  (
397  const UList<point>& points,
398  label& nTris,
399  label& nQuads
400  ) const;
401 
402  //- Split into triangles and quads.
403  // Results in triFaces (starting at triI) and quadFaces
404  // (starting at quadI).
405  // Returns number of new faces created.
406  label trianglesQuads
407  (
408  const UList<point>& points,
409  label& triI,
410  label& quadI,
411  faceList& triFaces,
412  faceList& quadFaces
413  ) const;
414 
415  //- Compare faces
416  // \return
417  // - 0: different
418  // - +1: identical
419  // - -1: same face, but different orientation
420  static int compare(const face& a, const face& b);
421 
422  //- Return true if the faces have the same vertices
423  static bool sameVertices(const face& a, const face& b);
424 };
425 
426 
427 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
428 
429 //- Hash specialization for face
430 template<>
431 struct Hash<face>
432 {
433  inline unsigned operator()(const face& obj, unsigned seed=0) const
434  {
435  return Hasher(obj.cdata(), obj.size()*sizeof(label), seed);
436  }
437 };
438 
439 
440 //- Specialization to offset faces, used in ListListOps::combineOffset
441 template<>
442 struct offsetOp<face>
443 {
444  inline face operator()
445  (
446  const face& x,
447  const label offset
448  ) const
449  {
450  face result(x.size());
451 
452  forAll(x, i)
453  {
454  result[i] = x[i] + offset;
455  }
456  return result;
457  }
458 };
459 
460 
461 //- Deprecated(2017-04) find the longest edge on a face.
462 //- Face point labels index into pts.
463 // \deprecated(2017-04) use class method instead
464 FOAM_DEPRECATED_FOR(2017-04, "use face::longestEdge() method")
465 label longestEdge(const face& f, const UList<point>& pts);
466 
467 
468 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
469 
470 inline bool operator==(const face& a, const face& b);
471 inline bool operator!=(const face& a, const face& b);
472 
473 
474 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
475 
476 } // End namespace Foam
477 
478 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
479 
480 #include "faceI.H"
481 
482 #ifdef NoRepository
483  #include "faceTemplates.C"
484 #endif
485 
486 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
487 
488 #endif
489 
490 // ************************************************************************* //
Foam::intersection::direction
direction
Definition: intersection.H:66
Foam::face::typeName
static const char *const typeName
Definition: face.H:143
Foam::Tensor< scalar >
Foam::face::reverseFace
face reverseFace() const
Return face with reverse direction.
Definition: face.C:605
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::face::points
pointField points(const UList< point > &points) const
Return the points corresponding to this face.
Definition: faceI.H:102
Foam::face::nearestPoint
pointHit nearestPoint(const point &p, const UList< point > &meshPoints) const
Return nearest point to face.
Definition: faceIntersection.C:200
intersection.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
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::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::Hasher
unsigned Hasher(const void *data, size_t len, unsigned seed=0)
Bob Jenkins's 96-bit mixer hashing function (lookup3)
Definition: Hasher.C:473
ListListOps.H
faceTemplates.C
Foam::face::mag
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition: faceI.H:130
Foam::face::nEdges
label nEdges() const
Return number of edges.
Definition: faceI.H:136
Foam::face::edgeDirection
int edgeDirection(const edge &e) const
The edge direction on the face.
Definition: face.C:760
Foam::face::faceEdge
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:143
Foam::face::inertia
tensor inertia(const UList< point > &p, const point &refPt=vector::zero, scalar density=1.0) const
Definition: face.C:707
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::intersection::FULL_RAY
Definition: intersection.H:74
Foam::face::NONE
Unknown proximity.
Definition: face.H:135
Foam::face::sameVertices
static bool sameVertices(const face &a, const face &b)
Return true if the faces have the same vertices.
Definition: face.C:402
Foam::face::longestEdge
label longestEdge(const UList< point > &pts) const
Find the longest edge on a face.
Definition: face.C:850
Foam::face::nTriangles
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:173
Foam::Hash::operator()
unsigned operator()(const T &obj, unsigned seed=0) const
Definition: Hash.H:57
pointHit.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::face::compare
static int compare(const face &a, const face &b)
Compare faces.
Definition: face.C:300
Foam::mode
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: MSwindows.C:564
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::face::unitNormal
vector unitNormal(const UList< point > &p) const
The unit normal.
Definition: faceI.H:122
Foam::face::trianglesQuads
label trianglesQuads(const UList< point > &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split into triangles and quads.
Definition: face.C:838
Foam::Hash
Hash function class. The default definition is for primitives, non-primitives used to hash entries on...
Definition: Hash.H:55
Foam::face::EDGE
Close to edge.
Definition: face.H:137
Foam::face::centre
point centre(const UList< point > &points) const
Centre point of face.
Definition: face.C:483
Foam::face::which
label which(const label pointLabel) const
Find local index on face for the point label,.
Definition: faceI.H:155
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
labelList.H
Foam::Field< vector >
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::Hash< face >::operator()
unsigned operator()(const face &obj, unsigned seed=0) const
Definition: face.H:432
Foam::intersection::algorithm
algorithm
Definition: intersection.H:72
Foam::face::ray
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
Definition: faceIntersection.C:36
Foam::face::average
Type average(const UList< point > &meshPoints, const Field< Type > &fld) const
Calculate average value at centroid of face.
Definition: faceTemplates.C:53
Foam::offsetOp
Offset operator for ListListOps::combineOffset()
Definition: ListListOps.H:101
Foam::face::proxType
proxType
The proximity classifications.
Definition: face.H:133
Foam::face::areaNormal
vector areaNormal(const UList< point > &p) const
The area normal - with magnitude equal to area of face.
Definition: face.C:547
fld
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))
Definition: gmvOutputLagrangian.H:23
Foam::face::intersection
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.
Definition: faceIntersection.C:142
edgeList.H
Foam::face::flip
void flip()
Flip the face in-place.
Definition: face.C:469
Foam::face::edges
edgeList edges() const
Return edges in face point ordering,.
Definition: face.C:742
Foam::face::areaInContact
scalar areaInContact(const UList< point > &meshPoints, const scalarField &v) const
Return area in contact, given the displacement in vertices.
Definition: faceAreaInContact.C:35
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::intersection::VECTOR
Definition: intersection.H:68
Foam::face::triangles
label triangles(const UList< point > &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition: face.C:810
pointField.H
Foam::face::prevLabel
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:167
Foam::face::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:161
Foam::face::nTrianglesQuads
label nTrianglesQuads(const UList< point > &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition: face.C:824
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
vectorField.H
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::UList
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:103
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::FOAM_DEPRECATED_FOR
class FOAM_DEPRECATED_FOR(2017-05, "Foam::Enum") NamedEnum
Definition: NamedEnum.H:69
Foam::face::POINT
Close to point.
Definition: face.H:136
Foam::face::FOAM_DEPRECATED_FOR
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal(const UList< point > &p) const
Legacy name for areaNormal()
Definition: face.H:213
Foam::face::contactSphereDiameter
scalar contactSphereDiameter(const point &p, const vector &n, const UList< point > &meshPoints) const
Return contact sphere diameter.
Definition: faceContactSphere.C:36
Foam::face::collapse
label collapse()
Collapse face by removing duplicate point labels.
Definition: face.C:444
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::zero
static const Vector< Cmpt > zero
Definition: VectorSpace.H:115
Foam::face::sign
int sign(const point &p, const UList< point > &points, const scalar tol=SMALL) const
Definition: faceIntersection.C:315
Foam::longestEdge
label longestEdge(const face &f, const UList< point > &pts)
Definition: face.C:872
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
faceListFwd.H
Forwards for various types of face lists.
Foam::face::face
constexpr face() noexcept=default
Default construct.
N
const Vector< label > N(dict.get< Vector< label >>("N"))
FixedList.H
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::face::nearestPointClassify
pointHit nearestPointClassify(const point &p, const UList< point > &meshPoints, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Definition: faceIntersection.C:214
triFace
face triFace(3)
Foam::face::found
bool found(const label pointLabel) const
Return true if the point label is found in face.
Definition: faceI.H:149
Foam::face::sweptVol
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:628