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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::face
29
30Description
31 A face is a list of labels corresponding to mesh vertices.
32
33See also
34 Foam::triFace
35
36SourceFiles
37 faceI.H
38 face.C
39 faceIntersection.C
40 faceContactSphere.C
41 faceAreaInContact.C
42 faceTemplates.C
43
44\*---------------------------------------------------------------------------*/
45
46#ifndef face_H
47#define face_H
48
49#include "pointField.H"
50#include "labelList.H"
51#include "edgeList.H"
52#include "vectorField.H"
53#include "faceListFwd.H"
54#include "intersection.H"
55#include "pointHit.H"
56#include "FixedList.H"
57#include "ListListOps.H"
58
59// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61namespace Foam
62{
63
64// Forward Declarations
65class face;
66class triFace;
67template<class T, int SizeMin> class DynamicList;
68
69/*---------------------------------------------------------------------------*\
70 Class face Declaration
71\*---------------------------------------------------------------------------*/
73class face
74:
75 public labelList
76{
77 // Private Member Functions
78
79 //- Construct list of edge vectors for face
80 tmp<vectorField> calcEdgeVectors
81 (
82 const UList<point>& points
83 ) const;
84
85 //- Find index of largest internal angle on face
86 label mostConcaveAngle
87 (
88 const UList<point>& points,
89 const vectorField& edges,
90 scalar& edgeCos
91 ) const;
92
93 //- Enumeration listing the modes for split()
94 enum splitMode
95 {
96 COUNTTRIANGLE,
97 COUNTQUAD,
98 SPLITTRIANGLE,
99 SPLITQUAD
100 };
101
102 //- Split face into triangles or triangles and quads.
103 // Stores results quadFaces[quadI], triFaces[triI]
104 // \return number of new faces created
105 label split
106 (
107 const splitMode mode,
108 const UList<point>& points,
109 label& triI,
110 label& quadI,
111 faceList& triFaces,
112 faceList& quadFaces
113 ) const;
114
115
116public:
117
118 // Public Typedefs
119
120 //- The proximity classifications
121 enum proxType
123 NONE = 0,
124 POINT,
126 };
127
128
129 // Static Data Members
131 static const char* const typeName;
132
133
134 // Constructors
135
136 //- Default construct
137 constexpr face() noexcept = default;
138
139 //- Construct given size, with invalid point labels (-1)
140 inline explicit face(const label sz);
141
142 //- Copy construct from list of point labels
143 inline explicit face(const labelUList& list);
144
145 //- Move construct from list of point labels
146 inline explicit face(labelList&& list);
147
148 //- Copy construct from an initializer list of point labels
149 inline explicit face(std::initializer_list<label> list);
150
151 //- Copy construct from list of point labels
152 template<unsigned N>
153 inline explicit face(const FixedList<label, N>& list);
154
155 //- Copy construct from subset of point labels
156 inline face(const labelUList& list, const labelUList& indices);
157
158 //- Copy construct from subset of point labels
159 template<unsigned N>
160 inline face(const labelUList& list, const FixedList<label, N>& indices);
161
162 //- Copy construct from triFace
163 face(const triFace& f);
164
165 //- Construct from Istream
166 inline face(Istream& is);
167
168
169 // Member Functions
170
171 //- Collapse face by removing duplicate point labels
172 // \return the collapsed size
173 label collapse();
174
175 //- Flip the face in-place.
176 // The starting points of the original and flipped face are identical.
177 void flip();
178
179 //- Return the points corresponding to this face
180 inline pointField points(const UList<point>& pts) const;
181
182 //- Centre point of face
183 point centre(const UList<point>& points) const;
184
185 //- Calculate average value at centroid of face
186 template<class Type>
187 Type average
188 (
189 const UList<point>& meshPoints,
190 const Field<Type>& fld
191 ) const;
192
193 //- The area normal - with magnitude equal to area of face
194 vector areaNormal(const UList<point>& p) const;
195
196 //- The unit normal
197 inline vector unitNormal(const UList<point>& p) const;
198
199 //- Legacy name for areaNormal()
200 // \deprecated(2018-06) Deprecated for new use
201 FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
202 vector normal(const UList<point>& p) const
203 {
204 return areaNormal(p); // Legacy definition
205 }
206
207 //- Magnitude of face area
208 inline scalar mag(const UList<point>& p) const;
209
210 //- Return face with reverse direction
211 // The starting points of the original and reverse face are identical.
212 face reverseFace() const;
213
214
215 // Navigation through face vertices
216
217 //- Find local index on face for the point label, same as find()
218 // \return position in face (0,1,2,...) or -1 if not found.
219 inline label which(const label pointLabel) const;
220
221 //- The vertex on face - identical to operator[], but with naming
222 //- similar to nextLabel(), prevLabel()
223 inline label thisLabel(const label i) const;
224
225 //- Next vertex on face
226 inline label nextLabel(const label i) const;
227
228 //- Previous vertex on face
229 inline label prevLabel(const label i) const;
230
231
232 //- Return the volume swept out by the face when its points move
233 scalar sweptVol
234 (
235 const UList<point>& oldPoints,
236 const UList<point>& newPoints
237 ) const;
238
239 //- Return the inertia tensor, with optional reference
240 //- point and density specification
242 (
243 const UList<point>& p,
244 const point& refPt = vector::zero,
245 scalar density = 1.0
246 ) const;
247
248 //- Return potential intersection with face with a ray starting
249 //- at p, direction n (does not need to be normalized)
250 // Does face-centre decomposition and returns triangle intersection
251 // point closest to p. Face-centre is calculated from point average.
252 // For a hit, the distance is signed. Positive number
253 // represents the point in front of triangle
254 // In case of miss the point is the nearest point on the face
255 // and the distance is the distance between the intersection point
256 // and the original point.
257 // The half-ray or full-ray intersection and the contact
258 // sphere adjustment of the projection vector is set by the
259 // intersection parameters
261 (
262 const point& p,
263 const vector& n,
264 const UList<point>& meshPoints,
267 ) const;
268
269 //- Fast intersection with a ray.
270 // Does face-centre decomposition and returns triangle intersection
271 // point closest to p. See triangle::intersection for details.
273 (
274 const point& p,
275 const vector& q,
276 const point& ctr,
277 const UList<point>& meshPoints,
278 const intersection::algorithm alg,
279 const scalar tol = 0.0
280 ) const;
281
282 //- Return nearest point to face
284 (
285 const point& p,
286 const UList<point>& meshPoints
287 ) const;
288
289 //- Return nearest point to face and classify it:
290 // + near point (nearType=POINT, nearLabel=0, 1, 2)
291 // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
292 // Note: edges are counted from starting vertex so
293 // e.g. edge n is from f[n] to f[0], where the face has n + 1
294 // points
296 (
297 const point& p,
298 const UList<point>& meshPoints,
299 label& nearType,
300 label& nearLabel
301 ) const;
302
303 //- The sign for the side of the face plane the point is on,
304 //- using three evenly distributed face points for the estimated normal.
305 // Uses the supplied tolerance for rounding around zero.
306 // \return
307 // - 0: on plane
308 // - +1: front-side
309 // - -1: back-side
310 int sign
311 (
312 const point& p,
313 const UList<point>& points,
314 const scalar tol = SMALL
315 ) const;
316
317 //- Return contact sphere diameter
319 (
320 const point& p,
321 const vector& n,
322 const UList<point>& meshPoints
323 ) const;
324
325 //- Return area in contact, given the displacement in vertices
326 scalar areaInContact
327 (
328 const UList<point>& meshPoints,
329 const scalarField& v
330 ) const;
331
332 //- Return number of edges
333 inline label nEdges() const noexcept;
334
335 //- Return i-th face edge (forward walk order).
336 // The edge 0 is from [0] to [1]
337 // and edge 1 is from [1] to [2]
338 inline Foam::edge edge(const label edgei) const;
339
340 //- Return vector of i-th face edge (forward walk order).
341 // The edge 0 is from [0] to [1]
342 // and edge 1 is from [1] to [2]
343 inline vector edge(const label edgei, const UList<point>& pts) const;
344
345 //- Return i-th face edge in reverse walk order.
346 // The rcEdge 0 is from [0] to [n-1]
347 // and rcEdge 1 is from [n-1] to [n-2]
348 inline Foam::edge rcEdge(const label edgei) const;
349
350 //- Return vector of i-th face edge in reverse walk order.
351 // The rcEdge 0 is from [0] to [n-1]
352 // and rcEdge 1 is from [n-1] to [n-2]
353 inline vector rcEdge(const label edgei, const UList<point>& pts) const;
354
355 //- Return list of edges in forward walk order.
356 // The edge 0 is from [0] to [1]
357 // and edge 1 is from [1] to [2]
358 edgeList edges() const;
359
360 //- Return list of edges in reverse walk order.
361 // The rcEdge 0 is from [0] to [n-1]
362 // and rcEdge 1 is from [n-1] to [n-2]
363 edgeList rcEdges() const;
364
365 //- Test the edge direction on the face
366 // \return
367 // - 0: edge not found on the face
368 // - +1: forward (counter-clockwise) on the face
369 // - -1: reverse (clockwise) on the face
370 int edgeDirection(const Foam::edge& e) const;
371
372 //- Find the longest edge on a face.
373 label longestEdge(const UList<point>& pts) const;
374
375
376 // Face splitting utilities
377
378 //- Number of triangles after splitting
379 inline label nTriangles() const;
380
381 //- Number of triangles after splitting
382 label nTriangles(const UList<point>& unused) const;
383
384 //- Split into triangles using existing points.
385 // Result in triFaces[triI..triI+nTri].
386 // Splits intelligently to maximize triangle quality.
387 // \Return number of faces created.
388 label triangles
389 (
390 const UList<point>& points,
391 label& triI,
392 faceList& triFaces
393 ) const;
394
395 //- Split into triangles using existing points.
396 // Append to DynamicList.
397 // \Return number of faces created.
398 template<int SizeMin>
399 label triangles
400 (
401 const UList<point>& points,
402 DynamicList<face, SizeMin>& triFaces
403 ) const;
404
405 //- Number of triangles and quads after splitting
406 // Returns the sum of both
407 label nTrianglesQuads
408 (
409 const UList<point>& points,
410 label& nTris,
411 label& nQuads
412 ) const;
413
414 //- Split into triangles and quads.
415 // Results in triFaces (starting at triI) and quadFaces
416 // (starting at quadI).
417 // Returns number of new faces created.
418 label trianglesQuads
419 (
420 const UList<point>& points,
421 label& triI,
422 label& quadI,
423 faceList& triFaces,
424 faceList& quadFaces
425 ) const;
426
427 //- Compare faces
428 // \return
429 // - 0: different
430 // - +1: identical
431 // - -1: same face, but different orientation
432 static int compare(const face& a, const face& b);
433
434 //- Return true if the faces have the same vertices
435 static bool sameVertices(const face& a, const face& b);
436
437
438 // Hashing
439
440 //- The symmetric hash value for face.
441 // Starts with lowest vertex value and walks in the direction
442 // of the next lowest value.
443 static unsigned symmhash_code(const UList<label>& f, unsigned seed=0);
444
445 //- The hash value for face.
446 // Currently hashes as sequential contiguous data,
447 // but could/should be commutative
448 inline unsigned hash_code(unsigned seed=0) const
449 {
450 return Foam::Hasher(this->cdata(), this->size_bytes(), seed);
451 }
452
453 //- The symmetric hash value for face.
454 // Starts with lowest vertex value and walks in the direction
455 // of the next lowest value.
456 inline unsigned symmhash_code(unsigned seed=0) const
457 {
458 return face::symmhash_code(*this, seed);
459 }
460
461 //- Hashing functor for face
462 struct hasher
464 unsigned operator()(const face& obj, unsigned seed=0) const
465 {
466 return obj.hash_code(seed);
467 }
468 };
469
470 //- Symmetric hashing functor for face
471 struct symmHasher
473 unsigned operator()(const face& obj, unsigned seed=0) const
474 {
475 return face::symmhash_code(obj, seed);
476 }
477 };
478
479
480 // Housekeeping
481
482 //- Identical to edge()
483 Foam::edge faceEdge(label edgei) const { return this->edge(edgei); }
484};
485
486
487// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
488
489//- Hash face as contiguous (non-commutative) list data
490template<> struct Hash<face> : face::hasher {};
491
492
493//- Specialization to offset faces, used in ListListOps::combineOffset
494template<>
495struct offsetOp<face>
497 inline face operator()
498 (
499 const face& x,
500 const label offset
501 ) const
502 {
503 face result(x.size());
504
505 forAll(x, i)
506 {
507 result[i] = x[i] + offset;
508 }
509 return result;
510 }
511};
512
513
514// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
515
516inline bool operator==(const face& a, const face& b);
517inline bool operator!=(const face& a, const face& b);
518
519
520// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
521
522} // End namespace Foam
523
524// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
525
526#include "faceI.H"
527
528#ifdef NoRepository
529 #include "faceTemplates.C"
530#endif
531
532// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
533
534#endif
535
536// ************************************************************************* //
label n
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:54
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
const T * cdata() const noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:230
std::streamsize size_bytes() const noexcept
Number of contiguous bytes for the List data.
Definition: UListI.H:258
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
label trianglesQuads(const UList< point > &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split into triangles and quads.
Definition: face.C:864
pointHit ray(const point &p, const vector &n, const UList< point > &meshPoints, const intersection::algorithm alg=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
pointHit nearestPointClassify(const point &p, const UList< point > &meshPoints, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
void flip()
Flip the face in-place.
Definition: face.C:500
static unsigned symmhash_code(const UList< label > &f, unsigned seed=0)
The symmetric hash value for face.
Definition: face.C:422
scalar sweptVol(const UList< point > &oldPoints, const UList< point > &newPoints) const
Return the volume swept out by the face when its points move.
Definition: face.C:659
label nEdges() const noexcept
Return number of edges.
Definition: faceI.H:118
label triangles(const UList< point > &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition: face.C:836
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition: faceI.H:112
unsigned hash_code(unsigned seed=0) const
The hash value for face.
Definition: face.H:447
label nTrianglesQuads(const UList< point > &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition: face.C:850
static int compare(const face &a, const face &b)
Compare faces.
Definition: face.C:281
tensor inertia(const UList< point > &p, const point &refPt=vector::zero, scalar density=1.0) const
Definition: face.C:738
scalar areaInContact(const UList< point > &meshPoints, const scalarField &v) const
Return area in contact, given the displacement in vertices.
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition: face.C:815
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition: faceI.H:141
label longestEdge(const UList< point > &pts) const
Find the longest edge on a face.
Definition: face.C:876
proxType
The proximity classifications.
Definition: face.H:121
@ POINT
Close to point.
Definition: face.H:123
@ EDGE
Close to edge.
Definition: face.H:124
@ NONE
Unknown proximity.
Definition: face.H:122
pointHit nearestPoint(const point &p, const UList< point > &meshPoints) const
Return nearest point to face.
point centre(const UList< point > &points) const
Centre point of face.
Definition: face.C:514
label which(const label pointLabel) const
Find local index on face for the point label, same as find()
Definition: faceI.H:163
scalar contactSphereDiameter(const point &p, const vector &n, const UList< point > &meshPoints) const
Return contact sphere diameter.
int sign(const point &p, const UList< point > &points, const scalar tol=SMALL) const
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:175
vector unitNormal(const UList< point > &p) const
The unit normal.
Definition: faceI.H:104
label thisLabel(const label i) const
Definition: faceI.H:169
edgeList rcEdges() const
Return list of edges in reverse walk order.
Definition: face.C:794
label collapse()
Collapse face by removing duplicate point labels.
Definition: face.C:475
face reverseFace() const
Return face with reverse direction.
Definition: face.C:636
Foam::edge faceEdge(label edgei) const
Identical to edge()
Definition: face.H:482
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:187
static bool sameVertices(const face &a, const face &b)
Return true if the faces have the same vertices.
Definition: face.C:382
edgeList edges() const
Return list of edges in forward walk order.
Definition: face.C:773
vector areaNormal(const UList< point > &p) const
The area normal - with magnitude equal to area of face.
Definition: face.C:578
Type average(const UList< point > &meshPoints, const Field< Type > &fld) const
Calculate average value at centroid of face.
Definition: faceTemplates.C:53
constexpr face() noexcept=default
Default construct.
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:181
static const char *const typeName
Definition: face.H:130
Foam::intersection.
Definition: intersection.H:53
A class for managing temporary objects.
Definition: tmp.H:65
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:72
volScalarField & p
Forwards for various types of face lists.
const pointField & points
Namespace for OpenFOAM.
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
unsigned Hasher(const void *data, size_t len, unsigned seed=0)
Bob Jenkins's 96-bit mixer hashing function (lookup3)
Definition: Hasher.C:581
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: MSwindows.C:572
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
const direction noexcept
Definition: Scalar.H:223
face triFace(3)
labelList f(nPoints)
volScalarField & b
Definition: createFields.H:27
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define FOAM_DEPRECATED_FOR(since, replacement)
Definition: stdFoam.H:52
Hash function class. The default definition is for primitives. Non-primitives used to hash entries on...
Definition: Hash.H:54
Hashing functor for face.
Definition: face.H:462
unsigned operator()(const face &obj, unsigned seed=0) const
Definition: face.H:463
Symmetric hashing functor for face.
Definition: face.H:471
unsigned operator()(const face &obj, unsigned seed=0) const
Definition: face.H:472
Offset operator for ListListOps::combineOffset()
Definition: ListListOps.H:102
A non-counting (dummy) refCount.
Definition: refCount.H:59
const Vector< label > N(dict.get< Vector< label > >("N"))