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-2023 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 "DynamicListFwd.H"
52 #include "intersection.H"
53 #include "pointHit.H"
54 #include "ListListOps.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 // Forward declaration of friend functions and operators
62 
63 class face;
64 class triFace;
65 
66 inline bool operator==(const face& a, const face& b);
67 inline bool operator!=(const face& a, const face& b);
68 inline Istream& operator>>(Istream&, face&);
69 
70 /*---------------------------------------------------------------------------*\
71  Class face Declaration
72 \*---------------------------------------------------------------------------*/
73 
74 class face
75 :
76  public labelList
77 {
78 public:
79 
80  //- Return types for classify
81  enum proxType
82  {
84  POINT, // Close to point
85  EDGE // Close to edge
86  };
87 
88  // Static Data Members
89 
90  static const char* const typeName;
91 
92 
93  // Static Functions
94 
95  //- Return vector area given face points
96  template<class PointField>
97  static vector area(const PointField& ps);
98 
99  //- Return centre point given face points
100  template<class PointField>
101  static vector centre(const PointField& ps);
102 
103  //- Return vector area and centre point given face points
104  template<class PointField>
106 
107  //- Return vector area and centre point given face points. Stabilised
108  // version. More expensive, but less likely to generate wildly
109  // inaccurate face-centres on degenerate faces.
110  template<class PointField>
112 
113 
114  // Constructors
115 
116  //- Construct null
117  inline face();
118 
119  //- Construct given size
120  explicit inline face(label);
121 
122  //- Construct from list of labels
123  explicit inline face(const labelUList&);
124 
125  //- Construct from list of labels
126  explicit inline face(const labelList&);
127 
128  //- Construct by transferring the parameter contents
129  explicit inline face(labelList&&);
130 
131  //- Copy construct from triFace
132  face(const triFace&);
133 
134  //- Construct from Istream
135  inline face(Istream&);
136 
137 
138  // Member Functions
139 
140  //- Collapse face by removing duplicate point labels
141  // return the collapsed size
142  label collapse();
143 
144  //- Flip the face in-place.
145  // The starting points of the original and flipped face are identical.
146  void flip();
147 
148  //- Return the points corresponding to this face
149  inline pointField points(const pointField&) const;
150 
151  //- Centre point of face
152  point centre(const pointField&) const;
153 
154  //- Calculate average value at centroid of face
155  template<class Type>
156  Type average(const pointField&, const Field<Type>&) const;
157 
158  //- Return vector area
159  vector area(const pointField&) const;
160 
161  //- Return scalar magnitude
162  inline scalar mag(const pointField&) const;
163 
164  //- Return unit normal
165  vector normal(const pointField&) const;
166 
167  //- Return face with reverse direction
168  // The starting points of the original and reverse face are identical.
169  face reverseFace() const;
170 
171  //- Navigation through face vertices
172 
173  //- Which vertex on face (face index given a global index)
174  // returns -1 if not found
175  label which(const label globalIndex) const;
176 
177  //- Next vertex on face
178  inline label nextLabel(const label i) const;
179 
180  //- Previous vertex on face
181  inline label prevLabel(const label i) const;
182 
183 
184  //- Return the volume swept out by the face when its points move
185  scalar sweptVol
186  (
187  const pointField& oldPoints,
188  const pointField& newPoints
189  ) const;
190 
191  //- Return the inertia tensor, with optional reference
192  // point and density specification
194  (
195  const pointField&,
196  const point& refPt = vector::zero,
197  scalar density = 1.0
198  ) const;
199 
200  //- Return potential intersection with face with a ray starting
201  // at p, direction n (does not need to be normalised)
202  // Does face-centre decomposition and returns triangle intersection
203  // point closest to p. Face-centre is calculated from point average.
204  // For a hit, the distance is signed. Positive number
205  // represents the point in front of triangle
206  // In case of miss the point is the nearest point on the face
207  // and the distance is the distance between the intersection point
208  // and the original point.
209  // The half-ray or full-ray intersection and the contact
210  // sphere adjustment of the projection vector is set by the
211  // intersection parameters
212  pointHit ray
213  (
214  const point& p,
215  const vector& n,
216  const pointField&,
217  const intersection::algorithm alg =
219  const intersection::direction dir =
221  ) const;
222 
223  //- Fast intersection with a ray.
224  // Does face-centre decomposition and returns triangle intersection
225  // point closest to p. See triangle::intersection for details.
227  (
228  const point& p,
229  const vector& q,
230  const point& ctr,
231  const pointField&,
232  const intersection::algorithm alg,
233  const scalar tol = 0.0
234  ) const;
235 
236  //- Return nearest point to face
238  (
239  const point& p,
240  const pointField&
241  ) const;
242 
243  //- Return nearest point to face and classify it:
244  // + near point (nearType=POINT, nearLabel=0, 1, 2)
245  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
246  // Note: edges are counted from starting vertex so
247  // e.g. edge n is from f[n] to f[0], where the face has n + 1
248  // points
250  (
251  const point& p,
252  const pointField&,
253  label& nearType,
254  label& nearLabel
255  ) const;
256 
257  //- Return contact sphere diameter
258  scalar contactSphereDiameter
259  (
260  const point& p,
261  const vector& n,
262  const pointField&
263  ) const;
264 
265  //- Return area in contact, given the displacement in vertices
266  scalar areaInContact
267  (
268  const pointField&,
269  const scalarField& v
270  ) const;
271 
272  //- Return number of edges
273  inline label nEdges() const;
274 
275  //- Return edges in face point ordering,
276  // i.e. edges()[0] is edge between [0] and [1]
277  edgeList edges() const;
278 
279  //- Return n-th face edge
280  inline edge faceEdge(const label n) const;
281 
282  //- Return the edge direction on the face
283  // Returns:
284  // - 0: edge not found on the face
285  // - +1: forward (counter-clockwise) on the face
286  // - -1: reverse (clockwise) on the face
287  int edgeDirection(const edge&) const;
288 
289  //- Size of the face's triangulation
290  inline label nTriangles() const;
291 
292  //- Compare faces
293  // 0: different
294  // +1: identical
295  // -1: same face, but different orientation
296  static int compare(const face&, const face&);
297 
298  //- Return true if the faces have the same vertices
299  static bool sameVertices(const face&, const face&);
300 
301 
302  // Member Operators
303 
304  //- Move assignment labelList
305  inline void operator=(labelList&&);
306 
307 
308  // Friend Operators
309 
310  friend bool operator==(const face& a, const face& b);
311  friend bool operator!=(const face& a, const face& b);
312 
313 
314  // Istream Operator
315 
316  friend Istream& operator>>(Istream&, face&);
317 };
318 
319 
320 //- Hash specialisation to offset faces in ListListOps::combineOffset
321 template<>
322 class offsetOp<face>
323 {
324 
325 public:
326 
327  inline face operator()
328  (
329  const face& x,
330  const label offset
331  ) const
332  {
333  face result(x.size());
334 
335  forAll(x, xI)
336  {
337  result[xI] = x[xI] + offset;
338  }
339  return result;
340  }
341 };
342 
343 
344 // Global functions
345 
346 //- Find the longest edge on a face. Face point labels index into pts.
347 label longestEdge(const face& f, const pointField& pts);
348 
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 } // End namespace Foam
353 
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355 
356 #include "faceI.H"
357 
358 #ifdef NoRepository
359  #include "faceTemplates.C"
360 #endif
361 
362 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
363 
364 #endif
365 
366 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
static const Form zero
Definition: VectorSpace.H:118
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
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.
void flip()
Flip the face in-place.
Definition: face.C:224
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: faceI.H:64
face()
Construct null.
Definition: faceI.H:28
tensor inertia(const pointField &, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: face.C:371
scalar areaInContact(const pointField &, const scalarField &v) const
Return area in contact, given the displacement in vertices.
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: faceI.H:81
int edgeDirection(const edge &) const
Return the edge direction on the face.
Definition: face.C:424
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.
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
friend Istream & operator>>(Istream &, face &)
static int compare(const face &, const face &)
Compare faces.
Definition: face.C:48
proxType
Return types for classify.
Definition: face.H:81
@ POINT
Definition: face.H:83
@ EDGE
Definition: face.H:84
@ NONE
Definition: face.H:82
static bool sameVertices(const face &, const face &)
Return true if the faces have the same vertices.
Definition: face.C:150
friend bool operator==(const face &a, const face &b)
scalar sweptVol(const pointField &oldPoints, const pointField &newPoints) const
Return the volume swept out by the face when its points move.
Definition: face.C:292
void operator=(labelList &&)
Move assignment labelList.
Definition: faceI.H:120
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:100
static vector centre(const PointField &ps)
Return centre point given face points.
vector normal(const pointField &) const
Return unit normal.
Definition: face.C:250
label collapse()
Collapse face by removing duplicate point labels.
Definition: face.C:199
face reverseFace() const
Return face with reverse direction.
Definition: face.C:256
label nTriangles() const
Size of the face's triangulation.
Definition: faceI.H:112
Type average(const pointField &, const Field< Type > &) const
Calculate average value at centroid of face.
edgeList edges() const
Return edges in face point ordering,.
Definition: face.C:406
static Tuple2< vector, point > areaAndCentreStabilised(const PointField &)
Return vector area and centre point given face points. Stabilised.
scalar contactSphereDiameter(const point &p, const vector &n, const pointField &) const
Return contact sphere diameter.
static vector area(const PointField &ps)
Return vector area given face points.
pointHit nearestPointClassify(const point &p, const pointField &, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
label which(const label globalIndex) const
Navigation through face vertices.
Definition: face.C:275
friend bool operator!=(const face &a, const face &b)
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:106
label nEdges() const
Return number of edges.
Definition: faceI.H:87
static const char *const typeName
Definition: face.H:89
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:94
static Tuple2< vector, point > areaAndCentre(const PointField &)
Return vector area and centre point given face points.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:71
volScalarField & b
Definition: createFields.H:25
Namespace for OpenFOAM.
bool operator!=(const particle &, const particle &)
Definition: particle.C:1210
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Istream & operator>>(Istream &, pistonPointEdgeData &)
label longestEdge(const face &f, const pointField &pts)
Find the longest edge on a face. Face point labels index into pts.
Definition: face.C:467
void offset(label &lst, const label o)
face triFace(3)
labelList f(nPoints)
volScalarField & p