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 #ifndef face_H
46 #define face_H
47 
48 #include "pointField.H"
49 #include "labelList.H"
50 #include "edgeList.H"
51 #include "vectorField.H"
52 #include "faceListFwd.H"
53 #include "intersection.H"
54 #include "pointHit.H"
55 #include "FixedList.H"
56 #include "ListListOps.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 
63 // Forward Declarations
64 class face;
65 class triFace;
66 
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  //- Return types for classify
132  enum proxType
133  {
135  POINT, // Close to point
136  EDGE // Close to edge
137  };
138 
139  // Static Data Members
140 
141  static const char* const typeName;
142 
143 
144  // Constructors
145 
146  //- Default construct
147  face() = default;
148 
149  //- Construct given size, with invalid point labels (-1)
150  inline explicit face(const label sz);
151 
152  //- Copy construct from list of labels
153  inline explicit face(const labelUList& list);
154 
155  //- Copy construct from list of labels
156  template<unsigned N>
157  inline explicit face(const FixedList<label, N>& list);
158 
159  //- Copy construct from an initializer list of labels
160  inline explicit face(std::initializer_list<label> list);
161 
162  //- Move construct from list of labels
163  inline explicit face(labelList&& list);
164 
165  //- Copy construct from triFace
166  face(const triFace& f);
167 
168  //- Construct from Istream
169  inline face(Istream& is);
170 
171 
172  // Member Functions
173 
174  //- Collapse face by removing duplicate point labels
175  // \return the collapsed size
176  label collapse();
177 
178  //- Flip the face in-place.
179  // The starting points of the original and flipped face are identical.
180  void flip();
181 
182  //- Return the points corresponding to this face
183  inline pointField points(const UList<point>& points) const;
184 
185  //- Centre point of face
186  point centre(const UList<point>& points) const;
187 
188  //- Calculate average value at centroid of face
189  template<class Type>
190  Type average
191  (
192  const UList<point>& meshPoints,
193  const Field<Type>& fld
194  ) const;
195 
196  //- The area normal - with magnitude equal to area of face
197  vector areaNormal(const UList<point>& p) const;
198 
199  //- The unit normal
200  inline vector unitNormal(const UList<point>& p) const;
201 
202  //- Legacy name for areaNormal()
203  // \deprecated(2018-06) Deprecated for new use
204  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
205  vector normal(const UList<point>& p) const
206  {
207  return areaNormal(p); // Legacy definition
208  }
209 
210  //- Magnitude of face area
211  inline scalar mag(const UList<point>& p) const;
212 
213  //- Return face with reverse direction
214  // The starting points of the original and reverse face are identical.
215  face reverseFace() const;
216 
217  // Navigation through face vertices
218 
219  //- Return true if the point label is found in face.
220  inline bool found(const label pointLabel) const;
221 
222  //- Find local index on face for the point label,
223  // \return position in face (0,1,2,...) or -1 if not found.
224  inline label which(const label pointLabel) const;
225 
226  //- Next vertex on face
227  inline label nextLabel(const label i) const;
228 
229  //- Previous vertex on face
230  inline label prevLabel(const label i) const;
231 
232 
233  //- Return the volume swept out by the face when its points move
234  scalar sweptVol
235  (
236  const UList<point>& oldPoints,
237  const UList<point>& newPoints
238  ) const;
239 
240  //- Return the inertia tensor, with optional reference
241  // point and density specification
243  (
244  const UList<point>& p,
245  const point& refPt = vector::zero,
246  scalar density = 1.0
247  ) const;
248 
249  //- Return potential intersection with face with a ray starting
250  // at p, direction n (does not need to be normalized)
251  // Does face-centre decomposition and returns triangle intersection
252  // point closest to p. Face-centre is calculated from point average.
253  // For a hit, the distance is signed. Positive number
254  // represents the point in front of triangle
255  // In case of miss the point is the nearest point on the face
256  // and the distance is the distance between the intersection point
257  // and the original point.
258  // The half-ray or full-ray intersection and the contact
259  // sphere adjustment of the projection vector is set by the
260  // intersection parameters
261  pointHit ray
262  (
263  const point& p,
264  const vector& n,
265  const UList<point>& meshPoints,
268  ) const;
269 
270  //- Fast intersection with a ray.
271  // Does face-centre decomposition and returns triangle intersection
272  // point closest to p. See triangle::intersection for details.
274  (
275  const point& p,
276  const vector& q,
277  const point& ctr,
278  const UList<point>& meshPoints,
279  const intersection::algorithm alg,
280  const scalar tol = 0.0
281  ) const;
282 
283  //- Return nearest point to face
285  (
286  const point& p,
287  const UList<point>& meshPoints
288  ) const;
289 
290  //- Return nearest point to face and classify it:
291  // + near point (nearType=POINT, nearLabel=0, 1, 2)
292  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
293  // Note: edges are counted from starting vertex so
294  // e.g. edge n is from f[n] to f[0], where the face has n + 1
295  // points
297  (
298  const point& p,
299  const UList<point>& meshPoints,
300  label& nearType,
301  label& nearLabel
302  ) const;
303 
304  //- The sign for the side of the face plane the point is on,
305  //- using three evenly distributed face points for the estimated normal.
306  // Uses the supplied tolerance for rounding around zero.
307  // \return
308  // - 0: on plane
309  // - +1: front-side
310  // - -1: back-side
311  int sign
312  (
313  const point& p,
314  const UList<point>& points,
315  const scalar tol = SMALL
316  ) const;
317 
318  //- Return contact sphere diameter
319  scalar contactSphereDiameter
320  (
321  const point& p,
322  const vector& n,
323  const UList<point>& meshPoints
324  ) const;
325 
326  //- Return area in contact, given the displacement in vertices
327  scalar areaInContact
328  (
329  const UList<point>& meshPoints,
330  const scalarField& v
331  ) const;
332 
333  //- Return number of edges
334  inline label nEdges() const;
335 
336  //- Return edges in face point ordering,
337  // i.e. edges()[0] is edge between [0] and [1]
338  edgeList edges() const;
339 
340  //- Return n-th face edge
341  inline edge faceEdge(const label n) const;
342 
343  //- The edge direction on the face
344  // \return
345  // - 0: edge not found on the face
346  // - +1: forward (counter-clockwise) on the face
347  // - -1: reverse (clockwise) on the face
348  int edgeDirection(const edge& e) const;
349 
350  //- Find the longest edge on a face.
351  label longestEdge(const UList<point>& pts) const;
352 
353  // Face splitting utilities
354 
355  //- Number of triangles after splitting
356  inline label nTriangles() const;
357 
358  //- Number of triangles after splitting
359  label nTriangles(const UList<point>& unused) const;
360 
361  //- Split into triangles using existing points.
362  // Result in triFaces[triI..triI+nTri].
363  // Splits intelligently to maximize triangle quality.
364  // \Return number of faces created.
365  label triangles
366  (
367  const UList<point>& points,
368  label& triI,
369  faceList& triFaces
370  ) const;
371 
372  //- Split into triangles using existing points.
373  // Append to DynamicList.
374  // \Return number of faces created.
375  template<int SizeMin>
376  label triangles
377  (
378  const UList<point>& points,
380  ) const;
381 
382  //- Number of triangles and quads after splitting
383  // Returns the sum of both
384  label nTrianglesQuads
385  (
386  const UList<point>& points,
387  label& nTris,
388  label& nQuads
389  ) const;
390 
391  //- Split into triangles and quads.
392  // Results in triFaces (starting at triI) and quadFaces
393  // (starting at quadI).
394  // Returns number of new faces created.
395  label trianglesQuads
396  (
397  const UList<point>& points,
398  label& triI,
399  label& quadI,
400  faceList& triFaces,
401  faceList& quadFaces
402  ) const;
403 
404  //- Compare faces
405  // \return
406  // - 0: different
407  // - +1: identical
408  // - -1: same face, but different orientation
409  static int compare(const face& a, const face& b);
410 
411  //- Return true if the faces have the same vertices
412  static bool sameVertices(const face& a, const face& b);
413 };
414 
415 
416 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
417 
418 //- Hash specialization for face
419 template<>
420 struct Hash<face>
421 {
422  inline unsigned operator()(const face& obj, unsigned seed=0) const
423  {
424  return Hasher(obj.cdata(), obj.size()*sizeof(label), seed);
425  }
426 };
427 
428 
429 //- Specialization to offset faces, used in ListListOps::combineOffset
430 template<>
431 struct offsetOp<face>
432 {
433  inline face operator()
434  (
435  const face& x,
436  const label offset
437  ) const
438  {
439  face result(x.size());
440 
441  forAll(x, i)
442  {
443  result[i] = x[i] + offset;
444  }
445  return result;
446  }
447 };
448 
449 
450 //- Deprecated(2017-04) find the longest edge on a face.
451 //- Face point labels index into pts.
452 // \deprecated(2017-04) use class method instead
453 label longestEdge(const face& f, const UList<point>& pts);
454 
455 
456 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
457 
458 inline bool operator==(const face& a, const face& b);
459 inline bool operator!=(const face& a, const face& b);
460 
461 
462 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
463 
464 } // End namespace Foam
465 
466 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
467 
468 #include "faceI.H"
469 
470 #ifdef NoRepository
471  #include "faceTemplates.C"
472 #endif
473 
474 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
475 
476 #endif
477 
478 // ************************************************************************* //
Foam::intersection::direction
direction
Definition: intersection.H:66
Foam::face::typeName
static const char *const typeName
Definition: face.H:140
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:85
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:59
Foam::PointHit< point >
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:113
Foam::face::nEdges
label nEdges() const
Return number of edges.
Definition: faceI.H:119
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:126
Foam::face::inertia
tensor inertia(const UList< point > &p, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
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
Definition: face.H:133
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:156
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
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:235
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::face::unitNormal
vector unitNormal(const UList< point > &p) const
The unit normal.
Definition: faceI.H:105
Foam::face::face
face()=default
Default construct.
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
Definition: face.H:135
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:138
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::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
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:421
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
Return potential intersection with face with a ray starting.
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
Return types for classify.
Definition: face.H:131
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:150
Foam::face::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:144
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:70
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
faceI.H
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::face::POINT
Definition: face.H:134
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:203
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< scalar >, scalar, 3 >::zero
static const Vector< scalar > 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.
FixedList.H
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:132
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