face.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::face
26 
27 Description
28  A face is a list of labels corresponding to mesh vertices.
29 
30 See also
31  Foam::triFace
32 
33 SourceFiles
34  faceI.H
35  face.C
36  faceIntersection.C
37  faceContactSphere.C
38  faceAreaInContact.C
39  faceTemplates.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef face_H
44 #define face_H
45 
46 #include "pointField.H"
47 #include "labelList.H"
48 #include "edgeList.H"
49 #include "vectorField.H"
50 #include "faceListFwd.H"
51 #include "intersection.H"
52 #include "pointHit.H"
53 #include "ListListOps.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 
60 // Forward declaration of friend functions and operators
61 
62 class face;
63 class triFace;
64 
65 template<class T, unsigned SizeInc, unsigned SizeMult, unsigned SizeDiv>
66 class DynamicList;
67 
68 inline bool operator==(const face& a, const face& b);
69 inline bool operator!=(const face& a, const face& b);
70 inline Istream& operator>>(Istream&, face&);
71 
72 /*---------------------------------------------------------------------------*\
73  Class face Declaration
74 \*---------------------------------------------------------------------------*/
75 
76 class face
77 :
78  public labelList
79 {
80  // Private Member Functions
81 
82  //- Edge to the right of face vertex i
83  inline label right(const label i) const;
84 
85  //- Edge to the left of face vertex i
86  inline label left(const label i) const;
87 
88  //- Construct list of edge vectors for face
89  tmp<vectorField> calcEdges
90  (
91  const pointField& points
92  ) const;
93 
94  //- Cos between neighbouring edges
95  scalar edgeCos
96  (
97  const vectorField& edges,
98  const label index
99  ) const;
100 
101  //- Find index of largest internal angle on face
102  label mostConcaveAngle
103  (
104  const pointField& points,
105  const vectorField& edges,
106  scalar& edgeCos
107  ) const;
108 
109  //- Enumeration listing the modes for split()
110  enum splitMode
111  {
112  COUNTTRIANGLE, // count if split into triangles
113  COUNTQUAD, // count if split into triangles&quads
114  SPLITTRIANGLE, // split into triangles
115  SPLITQUAD // split into triangles&quads
116  };
117 
118  //- Split face into triangles or triangles&quads.
119  // Stores results quadFaces[quadI], triFaces[triI]
120  // Returns number of new faces created
121  label split
122  (
123  const splitMode mode,
124  const pointField& points,
125  label& triI,
126  label& quadI,
127  faceList& triFaces,
128  faceList& quadFaces
129  ) const;
130 
131 
132 public:
133 
134  //- Return types for classify
135  enum proxType
136  {
138  POINT, // Close to point
139  EDGE // Close to edge
140  };
141 
142  // Static Data Members
144  static const char* const typeName;
145 
146 
147  // Static Functions
148 
149  //- Return centre point given face points
150  template<class PointField>
151  static vector centre(const PointField& ps);
152 
153  //- Return vector area given face points
154  template<class PointField>
155  static vector area(const PointField& ps);
156 
157 
158  // Constructors
159 
160  //- Construct null
161  inline face();
162 
163  //- Construct given size
164  explicit inline face(label);
165 
166  //- Construct from list of labels
167  explicit inline face(const labelUList&);
168 
169  //- Construct from list of labels
170  explicit inline face(const labelList&);
171 
172  //- Construct by transferring the parameter contents
173  explicit inline face(labelList&&);
174 
175  //- Copy construct from triFace
176  face(const triFace&);
177 
178  //- Construct from Istream
179  inline face(Istream&);
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 pointField&) const;
194 
195  //- Centre point of face
196  point centre(const pointField&) const;
197 
198  //- Calculate average value at centroid of face
199  template<class Type>
200  Type average(const pointField&, const Field<Type>&) const;
201 
202  //- Return vector area
203  vector area(const pointField&) const;
204 
205  //- Return scalar magnitude
206  inline scalar mag(const pointField&) const;
207 
208  //- Return unit normal
209  vector normal(const pointField&) const;
210 
211  //- Return face with reverse direction
212  // The starting points of the original and reverse face are identical.
213  face reverseFace() const;
214 
215  //- Navigation through face vertices
216 
217  //- Which vertex on face (face index given a global index)
218  // returns -1 if not found
219  label which(const label globalIndex) const;
220 
221  //- Next vertex on face
222  inline label nextLabel(const label i) const;
223 
224  //- Previous vertex on face
225  inline label prevLabel(const label i) const;
226 
227 
228  //- Return the volume swept out by the face when its points move
229  scalar sweptVol
230  (
231  const pointField& oldPoints,
232  const pointField& newPoints
233  ) const;
234 
235  //- Return the inertia tensor, with optional reference
236  // point and density specification
238  (
239  const pointField&,
240  const point& refPt = vector::zero,
241  scalar density = 1.0
242  ) const;
243 
244  //- Return potential intersection with face with a ray starting
245  // at p, direction n (does not need to be normalized)
246  // Does face-centre decomposition and returns triangle intersection
247  // point closest to p. Face-centre is calculated from point average.
248  // For a hit, the distance is signed. Positive number
249  // represents the point in front of triangle
250  // In case of miss the point is the nearest point on the face
251  // and the distance is the distance between the intersection point
252  // and the original point.
253  // The half-ray or full-ray intersection and the contact
254  // sphere adjustment of the projection vector is set by the
255  // intersection parameters
256  pointHit ray
257  (
258  const point& p,
259  const vector& n,
260  const pointField&,
261  const intersection::algorithm alg =
263  const intersection::direction dir =
265  ) const;
266 
267  //- Fast intersection with a ray.
268  // Does face-centre decomposition and returns triangle intersection
269  // point closest to p. See triangle::intersection for details.
271  (
272  const point& p,
273  const vector& q,
274  const point& ctr,
275  const pointField&,
276  const intersection::algorithm alg,
277  const scalar tol = 0.0
278  ) const;
279 
280  //- Return nearest point to face
282  (
283  const point& p,
284  const pointField&
285  ) const;
286 
287  //- Return nearest point to face and classify it:
288  // + near point (nearType=POINT, nearLabel=0, 1, 2)
289  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
290  // Note: edges are counted from starting vertex so
291  // e.g. edge n is from f[n] to f[0], where the face has n + 1
292  // points
294  (
295  const point& p,
296  const pointField&,
297  label& nearType,
298  label& nearLabel
299  ) const;
300 
301  //- Return contact sphere diameter
302  scalar contactSphereDiameter
303  (
304  const point& p,
305  const vector& n,
306  const pointField&
307  ) const;
308 
309  //- Return area in contact, given the displacement in vertices
310  scalar areaInContact
311  (
312  const pointField&,
313  const scalarField& v
314  ) const;
315 
316  //- Return number of edges
317  inline label nEdges() const;
318 
319  //- Return edges in face point ordering,
320  // i.e. edges()[0] is edge between [0] and [1]
321  edgeList edges() const;
322 
323  //- Return n-th face edge
324  inline edge faceEdge(const label n) const;
325 
326  //- Return the edge direction on the face
327  // Returns:
328  // - 0: edge not found on the face
329  // - +1: forward (counter-clockwise) on the face
330  // - -1: reverse (clockwise) on the face
331  int edgeDirection(const edge&) const;
332 
333  // Face splitting utilities
334 
335  //- Number of triangles after splitting
336  inline label nTriangles() const;
337 
338  //- Number of triangles after splitting
339  label nTriangles(const pointField& points) const;
340 
341  //- Split into triangles using existing points.
342  // Result in triFaces[triI..triI+nTri].
343  // Splits intelligently to maximize triangle quality.
344  // Returns number of faces created.
346  (
347  const pointField& points,
348  label& triI,
349  faceList& triFaces
350  ) const;
351 
352  //- Split into triangles using existing points.
353  // Append to DynamicList.
354  // Returns number of faces created.
355  template<unsigned SizeInc, unsigned SizeMult, unsigned SizeDiv>
357  (
358  const pointField& points,
360  ) const;
361 
362  //- Number of triangles and quads after splitting
363  // Returns the sum of both
365  (
366  const pointField& points,
367  label& nTris,
368  label& nQuads
369  ) const;
370 
371  //- Split into triangles and quads.
372  // Results in triFaces (starting at triI) and quadFaces
373  // (starting at quadI).
374  // Returns number of new faces created.
376  (
377  const pointField& points,
378  label& triI,
379  label& quadI,
380  faceList& triFaces,
381  faceList& quadFaces
382  ) const;
383 
384  //- Compare faces
385  // 0: different
386  // +1: identical
387  // -1: same face, but different orientation
388  static int compare(const face&, const face&);
389 
390  //- Return true if the faces have the same vertices
391  static bool sameVertices(const face&, const face&);
392 
393 
394  // Member Operators
395 
396  //- Move assignment labelList
397  inline void operator=(labelList&&);
398 
399 
400  // Friend Operators
401 
402  friend bool operator==(const face& a, const face& b);
403  friend bool operator!=(const face& a, const face& b);
404 
405 
406  // Istream Operator
407 
408  friend Istream& operator>>(Istream&, face&);
409 };
410 
411 
412 //- Hash specialization to offset faces in ListListOps::combineOffset
413 template<>
414 class offsetOp<face>
415 {
416 
417 public:
418 
419  inline face operator()
420  (
421  const face& x,
422  const label offset
423  ) const
424  {
425  face result(x.size());
426 
427  forAll(x, xI)
428  {
429  result[xI] = x[xI] + offset;
430  }
431  return result;
432  }
433 };
434 
435 
436 // Global functions
437 
438 //- Find the longest edge on a face. Face point labels index into pts.
439 label longestEdge(const face& f, const pointField& pts);
440 
441 
442 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
443 
444 } // End namespace Foam
445 
446 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
447 
448 #include "faceI.H"
449 
450 #ifdef NoRepository
451  #include "faceTemplates.C"
452 #endif
453 
454 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
455 
456 #endif
457 
458 // ************************************************************************* //
label nTrianglesQuads(const pointField &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition: face.C:738
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
mode_t mode(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file mode.
Definition: POSIX.C:461
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:59
friend Ostream & operator(Ostream &, const UList< T > &)
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
proxType
Return types for classify.
Definition: face.H:134
edgeList edges() const
Return edges in face point ordering,.
Definition: face.C:656
face()
Construct null.
Definition: faceI.H:44
friend Istream & operator>>(Istream &, face &)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label nEdges() const
Return number of edges.
Definition: faceI.H:103
pointHit ray(const point &p, const vector &n, const pointField &, const intersection::algorithm alg=intersection::algorithm::fullRay, const intersection::direction dir=intersection::direction::vector) const
Return potential intersection with face with a ray starting.
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:117
tensor inertia(const pointField &, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: face.C:621
static vector area(const PointField &ps)
Return vector area given face points.
label which(const label globalIndex) const
Navigation through face vertices.
Definition: face.C:525
label triangles(const pointField &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition: face.C:724
pointHit intersection(const point &p, const vector &q, const point &ctr, const pointField &, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:124
vector normal(const pointField &) const
Return unit normal.
Definition: face.C:500
label collapse()
Collapse face by removing duplicate point labels.
Definition: face.C:449
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
face reverseFace() const
Return face with reverse direction.
Definition: face.C:506
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
int edgeDirection(const edge &) const
Return the edge direction on the face.
Definition: face.C:674
const dimensionedScalar & b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
Istream & operator>>(Istream &, directionInfo &)
face triFace(3)
static int compare(const face &, const face &)
Compare faces.
Definition: face.C:298
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: faceI.H:97
label trianglesQuads(const pointField &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split into triangles and quads.
Definition: face.C:752
Type average(const pointField &, const Field< Type > &) const
Calculate average value at centroid of face.
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: faceI.H:80
scalar areaInContact(const pointField &, const scalarField &v) const
Return area in contact, given the displacement in vertices.
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:60
static const char *const typeName
Definition: face.H:143
static vector centre(const PointField &ps)
Return centre point given face points.
scalar sweptVol(const pointField &oldPoints, const pointField &newPoints) const
Return the volume swept out by the face when its points move.
Definition: face.C:542
label longestEdge(const face &f, const pointField &pts)
Find the longest edge on a face. Face point labels index into pts.
Definition: face.C:764
static bool sameVertices(const face &, const face &)
Return true if the faces have the same vertices.
Definition: face.C:400
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
void flip()
Flip the face in-place.
Definition: face.C:474
label n
friend bool operator==(const face &a, const face &b)
friend bool operator!=(const face &a, const face &b)
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:130
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:110
bool operator!=(const particle &, const particle &)
Definition: particle.C:1204
scalar contactSphereDiameter(const point &p, const vector &n, const pointField &) const
Return contact sphere diameter.
void operator=(labelList &&)
Move assignment labelList.
Definition: faceI.H:138
Namespace for OpenFOAM.
pointHit nearestPointClassify(const point &p, const pointField &, label &nearType, label &nearLabel) const
Return nearest point to face and classify it: