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 | Copyright (C) 2011-2015 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 SeeAlso
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
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
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
143 
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(const Xfer<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  //- Magnitude of face area
192  inline scalar mag(const pointField&) const;
193 
194  //- Vector normal; magnitude is equal to area of face
195  vector normal(const pointField&) const;
196 
197  //- Return face with reverse direction
198  // The starting points of the original and reverse face are identical.
199  face reverseFace() const;
200 
201  //- Navigation through face vertices
202 
203  //- Which vertex on face (face index given a global index)
204  // returns -1 if not found
205  label which(const label globalIndex) const;
206 
207  //- Next vertex on face
208  inline label nextLabel(const label i) const;
209 
210  //- Previous vertex on face
211  inline label prevLabel(const label i) const;
212 
213 
214  //- Return the volume swept out by the face when its points move
215  scalar sweptVol
216  (
217  const pointField& oldPoints,
218  const pointField& newPoints
219  ) const;
220 
221  //- Return the inertia tensor, with optional reference
222  // point and density specification
224  (
225  const pointField&,
226  const point& refPt = vector::zero,
227  scalar density = 1.0
228  ) const;
229 
230  //- Return potential intersection with face with a ray starting
231  // at p, direction n (does not need to be normalized)
232  // Does face-centre decomposition and returns triangle intersection
233  // point closest to p. Face-centre is calculated from point average.
234  // For a hit, the distance is signed. Positive number
235  // represents the point in front of triangle
236  // In case of miss the point is the nearest point on the face
237  // and the distance is the distance between the intersection point
238  // and the original point.
239  // The half-ray or full-ray intersection and the contact
240  // sphere adjustment of the projection vector is set by the
241  // intersection parameters
242  pointHit ray
243  (
244  const point& p,
245  const vector& n,
246  const pointField&,
249  ) const;
250 
251  //- Fast intersection with a ray.
252  // Does face-centre decomposition and returns triangle intersection
253  // point closest to p. See triangle::intersection for details.
255  (
256  const point& p,
257  const vector& q,
258  const point& ctr,
259  const pointField&,
260  const intersection::algorithm alg,
261  const scalar tol = 0.0
262  ) const;
263 
264  //- Return nearest point to face
266  (
267  const point& p,
268  const pointField&
269  ) const;
270 
271  //- Return nearest point to face and classify it:
272  // + near point (nearType=POINT, nearLabel=0, 1, 2)
273  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
274  // Note: edges are counted from starting vertex so
275  // e.g. edge n is from f[n] to f[0], where the face has n + 1
276  // points
278  (
279  const point& p,
280  const pointField&,
281  label& nearType,
282  label& nearLabel
283  ) const;
284 
285  //- Return contact sphere diameter
286  scalar contactSphereDiameter
287  (
288  const point& p,
289  const vector& n,
290  const pointField&
291  ) const;
292 
293  //- Return area in contact, given the displacement in vertices
294  scalar areaInContact
295  (
296  const pointField&,
297  const scalarField& v
298  ) const;
299 
300  //- Return number of edges
301  inline label nEdges() const;
302 
303  //- Return edges in face point ordering,
304  // i.e. edges()[0] is edge between [0] and [1]
305  edgeList edges() const;
306 
307  //- Return n-th face edge
308  inline edge faceEdge(const label n) const;
309 
310  //- Return the edge direction on the face
311  // Returns:
312  // - 0: edge not found on the face
313  // - +1: forward (counter-clockwise) on the face
314  // - -1: reverse (clockwise) on the face
315  int edgeDirection(const edge&) const;
316 
317  // Face splitting utilities
318 
319  //- Number of triangles after splitting
320  inline label nTriangles() const;
321 
322  //- Number of triangles after splitting
323  label nTriangles(const pointField& points) const;
324 
325  //- Split into triangles using existing points.
326  // Result in triFaces[triI..triI+nTri].
327  // Splits intelligently to maximize triangle quality.
328  // Returns number of faces created.
330  (
331  const pointField& points,
332  label& triI,
333  faceList& triFaces
334  ) const;
335 
336  //- Split into triangles using existing points.
337  // Append to DynamicList.
338  // Returns number of faces created.
339  template<unsigned SizeInc, unsigned SizeMult, unsigned SizeDiv>
341  (
342  const pointField& points,
344  ) const;
345 
346  //- Number of triangles and quads after splitting
347  // Returns the sum of both
349  (
350  const pointField& points,
351  label& nTris,
352  label& nQuads
353  ) const;
354 
355  //- Split into triangles and quads.
356  // Results in triFaces (starting at triI) and quadFaces
357  // (starting at quadI).
358  // Returns number of new faces created.
360  (
361  const pointField& points,
362  label& triI,
363  label& quadI,
364  faceList& triFaces,
365  faceList& quadFaces
366  ) const;
367 
368  //- Compare faces
369  // 0: different
370  // +1: identical
371  // -1: same face, but different orientation
372  static int compare(const face&, const face&);
373 
374  //- Return true if the faces have the same vertices
375  static bool sameVertices(const face&, const face&);
376 
377 
378  // Friend Operators
379 
380  friend bool operator==(const face& a, const face& b);
381  friend bool operator!=(const face& a, const face& b);
382 
383 
384  // Istream Operator
385 
386  friend Istream& operator>>(Istream&, face&);
387 };
388 
389 
390 //- Hash specialization to offset faces in ListListOps::combineOffset
391 template<>
392 class offsetOp<face>
393 {
394 
395 public:
396 
397  inline face operator()
398  (
399  const face& x,
400  const label offset
401  ) const
402  {
403  face result(x.size());
404 
405  forAll(x, xI)
406  {
407  result[xI] = x[xI] + offset;
408  }
409  return result;
410  }
411 };
412 
413 
414 // Global functions
415 
416 //- Find the longest edge on a face. Face point labels index into pts.
417 label longestEdge(const face& f, const pointField& pts);
418 
419 
420 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
421 
422 } // End namespace Foam
423 
424 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
425 
426 #include "faceI.H"
427 
428 #ifdef NoRepository
429 # include "faceTemplates.C"
430 #endif
431 
432 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
433 
434 #endif
435 
436 // ************************************************************************* //
Foam::longestEdge
label longestEdge(const face &f, const pointField &pts)
Find the longest edge on a face. Face point labels index into pts.
Definition: face.C:870
Foam::intersection::direction
direction
Definition: intersection.H:63
Foam::face::normal
vector normal(const pointField &) const
Vector normal; magnitude is equal to area of face.
Definition: face.C:552
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::face::typeName
static const char *const typeName
Definition: face.H:143
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::face::edgeDirection
int edgeDirection(const edge &) const
Return the edge direction on the face.
Definition: face.C:779
Foam::face::reverseFace
face reverseFace() const
Return face with reverse direction.
Definition: face.C:611
Foam::offsetOp
Definition: ListListOps.H:111
Foam::face::trianglesQuads
label trianglesQuads(const pointField &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split into triangles and quads.
Definition: face.C:858
Foam::face::ray
pointHit ray(const point &p, const vector &n, const pointField &, 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:43
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
p
p
Definition: pEqn.H:62
intersection.H
Foam::face::COUNTTRIANGLE
@ COUNTTRIANGLE
Definition: face.H:111
Foam::face::mag
scalar mag(const pointField &) const
Magnitude of face area.
Definition: faceI.H:97
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::face::average
Type average(const pointField &, const Field< Type > &) const
Calculate average value at centroid of face.
Definition: faceTemplates.C:51
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:51
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::face::operator!=
friend bool operator!=(const face &a, const face &b)
ListListOps.H
faceTemplates.C
Foam::face::nEdges
label nEdges() const
Return number of edges.
Definition: faceI.H:103
Foam::face::split
label split(const splitMode mode, const pointField &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split face into triangles or triangles&quads.
Definition: face.C:125
Foam::face::faceEdge
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:110
Foam::face::operator>>
friend Istream & operator>>(Istream &, face &)
Foam::face::centre
point centre(const pointField &) const
Centre point of face.
Definition: face.C:488
Foam::face::COUNTQUAD
@ COUNTQUAD
Definition: face.H:112
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::intersection::FULL_RAY
@ FULL_RAY
Definition: intersection.H:71
Foam::face::NONE
@ NONE
Definition: face.H:136
Foam::face::SPLITTRIANGLE
@ SPLITTRIANGLE
Definition: face.H:113
Foam::face::nTriangles
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:130
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::face::triangles
label triangles(const pointField &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition: face.C:830
Foam::face::contactSphereDiameter
scalar contactSphereDiameter(const point &p, const vector &n, const pointField &) const
Return contact sphere diameter.
Definition: faceContactSphere.C:34
Foam::face::edgeCos
scalar edgeCos(const vectorField &edges, const label index) const
Cos between neighbouring edges.
Definition: face.C:64
pointHit.H
Foam::Xfer
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::face::sweptVol
scalar sweptVol(const pointField &oldPoints, const pointField &newPoints) const
Return the volume swept out by the face when its points move.
Definition: face.C:647
Foam::mode
mode_t mode(const fileName &)
Return the file mode.
Definition: POSIX.C:573
Foam::face::which
label which(const label globalIndex) const
Navigation through face vertices.
Definition: face.C:630
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::face::areaInContact
scalar areaInContact(const pointField &, const scalarField &v) const
Return area in contact, given the displacement in vertices.
Definition: faceAreaInContact.C:36
Foam::face::EDGE
@ EDGE
Definition: face.H:138
Foam::face::calcEdges
tmp< vectorField > calcEdges(const pointField &points) const
Construct list of edge vectors for face.
Definition: face.C:41
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:59
Foam::face::nTrianglesQuads
label nTrianglesQuads(const pointField &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition: face.C:844
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
labelList.H
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::intersection::algorithm
algorithm
Definition: intersection.H:69
Foam::face::inertia
tensor inertia(const pointField &, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: face.C:726
Foam::face::mostConcaveAngle
label mostConcaveAngle(const pointField &points, const vectorField &edges, scalar &edgeCos) const
Find index of largest internal angle on face.
Definition: face.C:78
Foam::face::proxType
proxType
Return types for classify.
Definition: face.H:134
Foam::face::nearestPointClassify
pointHit nearestPointClassify(const point &p, const pointField &, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Definition: faceIntersection.C:213
Foam::operator!=
bool operator!=(const particle &, const particle &)
Definition: particle.C:147
edgeList.H
Foam::face::flip
void flip()
Flip the face in-place.
Definition: face.C:474
Foam::face::edges
edgeList edges() const
Return edges in face point ordering,.
Definition: face.C:761
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::intersection::VECTOR
@ VECTOR
Definition: intersection.H:65
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::face::operator==
friend bool operator==(const face &a, const face &b)
Foam::face::points
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: faceI.H:80
pointField.H
Foam::face::nearestPoint
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
Definition: faceIntersection.C:199
Foam::face::prevLabel
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:124
Foam::face::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:117
Foam::face::SPLITQUAD
@ SPLITQUAD
Definition: face.H:114
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
vectorField.H
Foam::face::right
label right(const label i) const
Edge to the right of face vertex i.
Definition: faceI.H:29
Foam::face::splitMode
splitMode
Enumeration listing the modes for split()
Definition: face.H:109
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:60
Foam::operator>>
Istream & operator>>(Istream &, edgeMesh &)
Definition: edgeMeshIO.C:141
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::face::POINT
@ POINT
Definition: face.H:137
Foam::face::collapse
label collapse()
Collapse face by removing duplicate point labels.
Definition: face.C:449
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::face::compare
static int compare(const face &, const face &)
Compare faces.
Definition: face.C:298
Foam::face::face
face()
Construct null.
Definition: faceI.H:44
Foam::face::left
label left(const label i) const
Edge to the left of face vertex i.
Definition: faceI.H:36
faceListFwd.H
Foam::face::sameVertices
static bool sameVertices(const face &, const face &)
Return true if the faces have the same vertices.
Definition: face.C:400
Foam::face::intersection
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.
Definition: faceIntersection.C:141
triFace
face triFace(3)