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-2019 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 inline Istream& operator>>(Istream& is, face& f);
70 
71 /*---------------------------------------------------------------------------*\
72  Class face Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 class face
76 :
77  public labelList
78 {
79  // Private Member Functions
80 
81  //- Edge to the right of face vertex i
82  inline label right(const label i) const;
83 
84  //- Edge to the left of face vertex i
85  inline label left(const label i) const;
86 
87  //- Construct list of edge vectors for face
88  tmp<vectorField> calcEdges
89  (
90  const UList<point>& points
91  ) const;
92 
93  //- Cos between neighbouring edges
94  scalar edgeCos
95  (
96  const vectorField& edges,
97  const label index
98  ) const;
99 
100  //- Find index of largest internal angle on face
101  label mostConcaveAngle
102  (
103  const UList<point>& points,
104  const vectorField& edges,
105  scalar& edgeCos
106  ) const;
107 
108  //- Enumeration listing the modes for split()
109  enum splitMode
110  {
111  COUNTTRIANGLE,
112  COUNTQUAD,
113  SPLITTRIANGLE,
114  SPLITQUAD
115  };
116 
117  //- Split face into triangles or triangles and quads.
118  // Stores results quadFaces[quadI], triFaces[triI]
119  // \return number of new faces created
120  label split
121  (
122  const splitMode mode,
123  const UList<point>& points,
124  label& triI,
125  label& quadI,
126  faceList& triFaces,
127  faceList& quadFaces
128  ) const;
129 
130 
131 public:
132 
133  //- Return types for classify
134  enum proxType
135  {
137  POINT, // Close to point
138  EDGE // Close to edge
139  };
140 
141  // Static data members
142 
143  static const char* const typeName;
144 
145 
146  // Constructors
147 
148  //- Construct null
149  inline face();
150 
151  //- Construct given size, with invalid point labels (-1)
152  inline explicit face(const label sz);
153 
154  //- Copy construct from list of labels
155  inline explicit face(const labelUList& list);
156 
157  //- Copy construct from list of labels
158  template<unsigned N>
159  inline explicit face(const FixedList<label, N>& list);
160 
161  //- Copy construct from an initializer list of labels
162  inline explicit face(std::initializer_list<label> list);
163 
164  //- Move construct from list of labels
165  inline explicit face(labelList&& list);
166 
167  //- Copy construct from triFace
168  face(const triFace& f);
169 
170  //- Construct from Istream
171  inline face(Istream& is);
172 
173 
174  // Member Functions
175 
176  //- Collapse face by removing duplicate point labels
177  // return the collapsed size
178  label collapse();
179 
180  //- Flip the face in-place.
181  // The starting points of the original and flipped face are identical.
182  void flip();
183 
184  //- Return the points corresponding to this face
185  inline pointField points(const UList<point>& points) const;
186 
187  //- Centre point of face
188  point centre(const UList<point>& points) const;
189 
190  //- Calculate average value at centroid of face
191  template<class Type>
192  Type average
193  (
194  const UList<point>& meshPoints,
195  const Field<Type>& fld
196  ) const;
197 
198  //- The area normal - with magnitude equal to area of face
199  vector areaNormal(const UList<point>& p) const;
200 
201  //- The unit normal
202  inline vector unitNormal(const UList<point>& p) const;
203 
204  //- Legacy name for areaNormal()
205  // \deprecated(2018-06) Deprecated for new use
206  vector FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
207  normal(const UList<point>& p) const
208  {
209  return areaNormal(p); // Legacy definition
210  }
211 
212  //- Magnitude of face area
213  inline scalar mag(const UList<point>& p) const;
214 
215  //- Return face with reverse direction
216  // The starting points of the original and reverse face are identical.
217  face reverseFace() const;
218 
219  // Navigation through face vertices
220 
221  //- Return true if the point label is found in face.
222  inline bool found(const label pointLabel) const;
223 
224  //- Find local index on face for the point label,
225  // \return position in face (0,1,2,...) or -1 if not found.
226  inline label which(const label pointLabel) 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;
337 
338  //- Return edges in face point ordering,
339  // i.e. edges()[0] is edge between [0] and [1]
340  edgeList edges() const;
341 
342  //- Return n-th face edge
343  inline edge faceEdge(const label n) const;
344 
345  //- The edge direction on the face
346  // \return
347  // - 0: edge not found on the face
348  // - +1: forward (counter-clockwise) on the face
349  // - -1: reverse (clockwise) on the face
350  int edgeDirection(const edge& e) const;
351 
352  //- Find the longest edge on a face.
353  label longestEdge(const UList<point>& pts) const;
354 
355  // Face splitting utilities
356 
357  //- Number of triangles after splitting
358  inline label nTriangles() const;
359 
360  //- Number of triangles after splitting
361  label nTriangles(const UList<point>& unused) const;
362 
363  //- Split into triangles using existing points.
364  // Result in triFaces[triI..triI+nTri].
365  // Splits intelligently to maximize triangle quality.
366  // Returns number of faces created.
368  (
369  const UList<point>& points,
370  label& triI,
371  faceList& triFaces
372  ) const;
373 
374  //- Split into triangles using existing points.
375  // Append to DynamicList.
376  // Returns number of faces created.
377  template<int SizeMin>
379  (
380  const UList<point>& points,
382  ) const;
383 
384  //- Number of triangles and quads after splitting
385  // Returns the sum of both
387  (
388  const UList<point>& points,
389  label& nTris,
390  label& nQuads
391  ) const;
392 
393  //- Split into triangles and quads.
394  // Results in triFaces (starting at triI) and quadFaces
395  // (starting at quadI).
396  // Returns number of new faces created.
398  (
399  const UList<point>& points,
400  label& triI,
401  label& quadI,
402  faceList& triFaces,
403  faceList& quadFaces
404  ) const;
405 
406  //- Compare faces
407  // \return
408  // - 0: different
409  // - +1: identical
410  // - -1: same face, but different orientation
411  static int compare(const face& a, const face& b);
412 
413  //- Return true if the faces have the same vertices
414  static bool sameVertices(const face& a, const face& b);
415 
416 
417  // Istream Operator
418 
419  friend Istream& operator>>(Istream& is, face& f);
420 };
421 
422 
423 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
424 
425 //- Hash specialization for face
426 template<>
427 struct Hash<face>
428 {
429  inline unsigned operator()(const face& obj, unsigned seed=0) const
430  {
431  return Hasher(obj.cdata(), obj.size()*sizeof(label), seed);
432  }
433 };
434 
435 
436 //- Specialization to offset faces, used in ListListOps::combineOffset
437 template<>
438 struct offsetOp<face>
439 {
440  inline face operator()
441  (
442  const face& x,
443  const label offset
444  ) const
445  {
446  face result(x.size());
447 
448  forAll(x, i)
449  {
450  result[i] = x[i] + offset;
451  }
452  return result;
453  }
454 };
455 
456 
457 //- Deprecated(2017-04) find the longest edge on a face.
458 //- Face point labels index into pts.
459 // \deprecated(2017-04) use class method instead
460 label longestEdge(const face& f, const UList<point>& pts);
461 
462 
463 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
464 
465 inline bool operator==(const face& a, const face& b);
466 inline bool operator!=(const face& a, const face& b);
467 
468 
469 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
470 
471 } // End namespace Foam
472 
473 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
474 
475 #include "faceI.H"
476 
477 #ifdef NoRepository
478  #include "faceTemplates.C"
479 #endif
480 
481 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
482 
483 #endif
484 
485 // ************************************************************************* //
Foam::intersection::direction
direction
Definition: intersection.H:66
Foam::face::typeName
static const char *const typeName
Definition: face.H:142
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:91
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:57
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:119
Foam::face::nEdges
label nEdges() const
Return number of edges.
Definition: faceI.H:125
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:132
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: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:162
Foam::Hash::operator()
unsigned operator()(const T &obj, unsigned seed=0) const
Definition: Hash.H:57
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:228
pointHit.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
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:111
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:137
Foam::face::operator>>
friend Istream & operator>>(Istream &is, face &f)
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
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:144
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:428
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: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
Foam::face::FOAM_DEPRECATED_FOR
vector FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") normal(const UList< point > &p) const
Legacy name for areaNormal()
Definition: face.H:205
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:156
Foam::face::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:150
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:136
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:74
Foam::face::face
face()
Construct null.
Definition: faceI.H:47
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:138
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