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-2019 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  // Constructors
148 
149  //- Construct null
150  inline face();
151 
152  //- Construct given size
153  explicit inline face(label);
154 
155  //- Construct from list of labels
156  explicit inline face(const labelUList&);
157 
158  //- Construct from list of labels
159  explicit inline face(const labelList&);
160 
161  //- Construct by transferring the parameter contents
162  explicit inline face(labelList&&);
163 
164  //- Copy construct from triFace
165  face(const triFace&);
166 
167  //- Construct from Istream
168  inline face(Istream&);
169 
170 
171  // Member Functions
172 
173  //- Collapse face by removing duplicate point labels
174  // return the collapsed size
175  label collapse();
176 
177  //- Flip the face in-place.
178  // The starting points of the original and flipped face are identical.
179  void flip();
180 
181  //- Return the points corresponding to this face
182  inline pointField points(const pointField&) const;
183 
184  //- Centre point of face
185  point centre(const pointField&) const;
186 
187  //- Calculate average value at centroid of face
188  template<class Type>
189  Type average(const pointField&, const Field<Type>&) const;
190 
191  //- Return vector area
192  vector area(const pointField&) const;
193 
194  //- Return scalar magnitude
195  inline scalar mag(const pointField&) const;
196 
197  //- Return unit normal
198  vector normal(const pointField&) const;
199 
200  //- Return face with reverse direction
201  // The starting points of the original and reverse face are identical.
202  face reverseFace() const;
203 
204  //- Navigation through face vertices
205 
206  //- Which vertex on face (face index given a global index)
207  // returns -1 if not found
208  label which(const label globalIndex) const;
209 
210  //- Next vertex on face
211  inline label nextLabel(const label i) const;
212 
213  //- Previous vertex on face
214  inline label prevLabel(const label i) const;
215 
216 
217  //- Return the volume swept out by the face when its points move
218  scalar sweptVol
219  (
220  const pointField& oldPoints,
221  const pointField& newPoints
222  ) const;
223 
224  //- Return the inertia tensor, with optional reference
225  // point and density specification
227  (
228  const pointField&,
229  const point& refPt = vector::zero,
230  scalar density = 1.0
231  ) const;
232 
233  //- Return potential intersection with face with a ray starting
234  // at p, direction n (does not need to be normalized)
235  // Does face-centre decomposition and returns triangle intersection
236  // point closest to p. Face-centre is calculated from point average.
237  // For a hit, the distance is signed. Positive number
238  // represents the point in front of triangle
239  // In case of miss the point is the nearest point on the face
240  // and the distance is the distance between the intersection point
241  // and the original point.
242  // The half-ray or full-ray intersection and the contact
243  // sphere adjustment of the projection vector is set by the
244  // intersection parameters
245  pointHit ray
246  (
247  const point& p,
248  const vector& n,
249  const pointField&,
250  const intersection::algorithm alg =
252  const intersection::direction dir =
254  ) const;
255 
256  //- Fast intersection with a ray.
257  // Does face-centre decomposition and returns triangle intersection
258  // point closest to p. See triangle::intersection for details.
260  (
261  const point& p,
262  const vector& q,
263  const point& ctr,
264  const pointField&,
265  const intersection::algorithm alg,
266  const scalar tol = 0.0
267  ) const;
268 
269  //- Return nearest point to face
271  (
272  const point& p,
273  const pointField&
274  ) const;
275 
276  //- Return nearest point to face and classify it:
277  // + near point (nearType=POINT, nearLabel=0, 1, 2)
278  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
279  // Note: edges are counted from starting vertex so
280  // e.g. edge n is from f[n] to f[0], where the face has n + 1
281  // points
283  (
284  const point& p,
285  const pointField&,
286  label& nearType,
287  label& nearLabel
288  ) const;
289 
290  //- Return contact sphere diameter
291  scalar contactSphereDiameter
292  (
293  const point& p,
294  const vector& n,
295  const pointField&
296  ) const;
297 
298  //- Return area in contact, given the displacement in vertices
299  scalar areaInContact
300  (
301  const pointField&,
302  const scalarField& v
303  ) const;
304 
305  //- Return number of edges
306  inline label nEdges() const;
307 
308  //- Return edges in face point ordering,
309  // i.e. edges()[0] is edge between [0] and [1]
310  edgeList edges() const;
311 
312  //- Return n-th face edge
313  inline edge faceEdge(const label n) const;
314 
315  //- Return the edge direction on the face
316  // Returns:
317  // - 0: edge not found on the face
318  // - +1: forward (counter-clockwise) on the face
319  // - -1: reverse (clockwise) on the face
320  int edgeDirection(const edge&) const;
321 
322  // Face splitting utilities
323 
324  //- Number of triangles after splitting
325  inline label nTriangles() const;
326 
327  //- Number of triangles after splitting
328  label nTriangles(const pointField& points) const;
329 
330  //- Split into triangles using existing points.
331  // Result in triFaces[triI..triI+nTri].
332  // Splits intelligently to maximize triangle quality.
333  // Returns number of faces created.
335  (
336  const pointField& points,
337  label& triI,
338  faceList& triFaces
339  ) const;
340 
341  //- Split into triangles using existing points.
342  // Append to DynamicList.
343  // Returns number of faces created.
344  template<unsigned SizeInc, unsigned SizeMult, unsigned SizeDiv>
346  (
347  const pointField& points,
349  ) const;
350 
351  //- Number of triangles and quads after splitting
352  // Returns the sum of both
354  (
355  const pointField& points,
356  label& nTris,
357  label& nQuads
358  ) const;
359 
360  //- Split into triangles and quads.
361  // Results in triFaces (starting at triI) and quadFaces
362  // (starting at quadI).
363  // Returns number of new faces created.
365  (
366  const pointField& points,
367  label& triI,
368  label& quadI,
369  faceList& triFaces,
370  faceList& quadFaces
371  ) const;
372 
373  //- Compare faces
374  // 0: different
375  // +1: identical
376  // -1: same face, but different orientation
377  static int compare(const face&, const face&);
378 
379  //- Return true if the faces have the same vertices
380  static bool sameVertices(const face&, const face&);
381 
382 
383  // Member Operators
384 
385  //- Move assignment labelList
386  inline void operator=(labelList&&);
387 
388 
389  // Friend Operators
390 
391  friend bool operator==(const face& a, const face& b);
392  friend bool operator!=(const face& a, const face& b);
393 
394 
395  // Istream Operator
396 
397  friend Istream& operator>>(Istream&, face&);
398 };
399 
400 
401 //- Hash specialization to offset faces in ListListOps::combineOffset
402 template<>
403 class offsetOp<face>
404 {
405 
406 public:
407 
408  inline face operator()
409  (
410  const face& x,
411  const label offset
412  ) const
413  {
414  face result(x.size());
415 
416  forAll(x, xI)
417  {
418  result[xI] = x[xI] + offset;
419  }
420  return result;
421  }
422 };
423 
424 
425 // Global functions
426 
427 //- Find the longest edge on a face. Face point labels index into pts.
428 label longestEdge(const face& f, const pointField& pts);
429 
430 
431 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
432 
433 } // End namespace Foam
434 
435 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
436 
437 #include "faceI.H"
438 
439 #ifdef NoRepository
440  #include "faceTemplates.C"
441 #endif
442 
443 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
444 
445 #endif
446 
447 // ************************************************************************* //
label nTrianglesQuads(const pointField &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition: face.C:851
#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:769
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:734
point centre(const pointField &) const
Centre point of face.
Definition: face.C:488
label which(const label globalIndex) const
Navigation through face vertices.
Definition: face.C:638
label triangles(const pointField &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition: face.C:837
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:611
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:619
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:787
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:865
Type average(const pointField &, const Field< Type > &) const
Calculate average value at centroid of face.
Definition: faceTemplates.C:51
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
vector area(const pointField &) const
Return vector area.
Definition: face.C:552
static const char *const typeName
Definition: face.H:143
scalar sweptVol(const pointField &oldPoints, const pointField &newPoints) const
Return the volume swept out by the face when its points move.
Definition: face.C:655
label longestEdge(const face &f, const pointField &pts)
Find the longest edge on a face. Face point labels index into pts.
Definition: face.C:877
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:1223
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: