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