Loading...
Searching...
No Matches
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-2025 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 Foam_face_H
47#define Foam_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;
69/*---------------------------------------------------------------------------*\
70 Class face Declaration
71\*---------------------------------------------------------------------------*/
72
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
122 {
123 NONE = 0,
124 POINT,
125 EDGE
126 };
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 vertex labels (-1)
140 inline explicit face(const label sz);
141
142 //- Copy construct from list of vertex labels
143 inline explicit face(const labelUList& list);
144
145 //- Move construct from list of vertex labels
146 inline explicit face(labelList&& list);
147
148 //- Copy construct from an initializer list of vertex labels
149 inline explicit face(std::initializer_list<label> list);
150
151 //- Copy construct from list of vertex labels
152 template<unsigned N>
153 inline explicit face(const FixedList<label, N>& list);
154
155 //- Copy construct from subset of vertex labels
156 inline face(const labelUList& list, const labelUList& indices);
157
158 //- Copy construct from subset of vertex 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 vertex labels
172 // \return the collapsed size
173 label collapse();
174
175 //- Flip the face in-place.
176 // The starting vertex 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 //- Average of face points: a quick estimate of the face centre.
186 inline point average(const UList<point>& points) const;
187
188 //- Calculate average value at centroid of face
189 template<class Type>
190 Type average
191 (
192 const UList<point>& meshPoints,
193 const Field<Type>& fld
194 ) const;
195
196 //- The area normal - with magnitude equal to area of face
197 vector areaNormal(const UList<point>& p) const;
198
199 //- The unit normal
200 inline vector unitNormal(const UList<point>& p) const;
201
202 //- Legacy name for areaNormal()
203 // \deprecated(2018-06) Deprecated for new use
204 FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
205 vector normal(const UList<point>& p) const
206 {
207 return areaNormal(p); // Legacy definition
208 }
209
210 //- Magnitude of face area
211 inline scalar mag(const UList<point>& p) const;
212
213 //- Magnitude squared of face area
214 inline scalar magSqr(const UList<point>& p) const;
215
216 //- The enclosing (bounding) box for the face
217 inline Pair<point> box(const UList<point>& pts) const;
218
219 //- Return face with reverse direction
220 // The starting vertex of the original and reverse face are identical.
221 face reverseFace() const;
222
223
224 // Navigation through face vertices
225
226 //- Find local vertex on face for the vertex label, same as find()
227 // \return position in face (0,1,2,...) or -1 if not found.
228 inline label which(const label vertex) const;
229
230 //- The vertex on face - identical to operator[], but with naming
231 //- similar to nextLabel(), prevLabel()
232 inline label thisLabel(const label i) const;
233
234 //- Next vertex on face
235 inline label nextLabel(const label i) const;
236
237 //- Previous vertex on face
238 inline label prevLabel(const label i) const;
239
240
241 //- Return the volume swept out by the face when its points move
242 scalar sweptVol
243 (
244 const UList<point>& oldPoints,
245 const UList<point>& newPoints
246 ) const;
247
248 //- Return the inertia tensor, with optional reference
249 //- point and density specification
252 const UList<point>& p,
253 const point& refPt = vector::zero,
254 scalar density = 1.0
255 ) const;
256
257 //- Return potential intersection with face with a ray starting
258 //- at p, direction n (does not need to be normalized)
259 // Does face-centre decomposition and returns triangle intersection
260 // point closest to p. Face-centre is calculated from point average.
261 // For a hit, the distance is signed. Positive number
262 // represents the point in front of triangle
263 // In case of miss the point is the nearest point on the face
264 // and the distance is the distance between the intersection point
265 // and the original point.
266 // The half-ray or full-ray intersection and the contact
267 // sphere adjustment of the projection vector is set by the
268 // intersection parameters
270 (
271 const point& p,
272 const vector& n,
273 const UList<point>& meshPoints,
276 ) const;
277
278 //- Fast intersection with a ray.
279 // Does face-centre decomposition and returns triangle intersection
280 // point closest to p. See triangle::intersection for details.
282 (
283 const point& p,
284 const vector& q,
285 const point& ctr,
286 const UList<point>& meshPoints,
287 const intersection::algorithm alg,
288 const scalar tol = 0.0
289 ) const;
290
291 //- Return nearest point to face
293 (
294 const point& p,
295 const UList<point>& meshPoints
296 ) const;
297
298 //- Return nearest point to face and classify it:
299 // + near point (nearType=POINT, nearLabel=0, 1, 2)
300 // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
301 // Note: edges are counted from starting vertex so
302 // e.g. edge n is from f[n] to f[0], where the face has n + 1
303 // points
305 (
306 const point& p,
307 const UList<point>& meshPoints,
308 label& nearType,
309 label& nearLabel
310 ) const;
311
312 //- The sign for the side of the face plane the point is on,
313 //- using three evenly distributed face points for the estimated normal.
314 // Uses the supplied tolerance for rounding around zero.
315 // \return
316 // - 0: on plane
317 // - +1: front-side
318 // - -1: back-side
319 int sign
320 (
321 const point& p,
322 const UList<point>& points,
323 const scalar tol = SMALL
324 ) const;
325
326 //- Return contact sphere diameter
328 (
329 const point& p,
330 const vector& n,
331 const UList<point>& meshPoints
332 ) const;
333
334 //- Return area in contact, given the displacement in vertices
335 scalar areaInContact
336 (
337 const UList<point>& meshPoints,
338 const scalarField& v
339 ) const;
340
341 //- Return number of edges
342 inline label nEdges() const noexcept;
343
344 //- Return i-th face edge (forward walk order).
345 // The edge 0 is from [0] to [1]
346 // and edge 1 is from [1] to [2]
347 inline Foam::edge edge(const label edgei) const;
348
349 //- Return vector of i-th face edge (forward walk order).
350 // The edge 0 is from [0] to [1]
351 // and edge 1 is from [1] to [2]
352 inline vector edge(const label edgei, const UList<point>& pts) const;
353
354 //- Return i-th face edge in reverse walk order.
355 // The rcEdge 0 is from [0] to [n-1]
356 // and rcEdge 1 is from [n-1] to [n-2]
357 inline Foam::edge rcEdge(const label edgei) const;
358
359 //- Return vector of i-th face edge in reverse walk order.
360 // The rcEdge 0 is from [0] to [n-1]
361 // and rcEdge 1 is from [n-1] to [n-2]
362 inline vector rcEdge(const label edgei, const UList<point>& pts) const;
363
364 //- Return list of edges in forward walk order.
365 // The edge 0 is from [0] to [1]
366 // and edge 1 is from [1] to [2]
367 edgeList edges() const;
368
369 //- Return list of edges in reverse walk order.
370 // The rcEdge 0 is from [0] to [n-1]
371 // and rcEdge 1 is from [n-1] to [n-2]
372 edgeList rcEdges() const;
373
374
375 // Searching
376
377 //- Regular contains(pointi) tests
378 using labelList::contains;
379
380 //- True if face contains(edge)
381 inline bool contains(const Foam::edge& e) const;
382
383 //- Regular find pointi within face
384 using labelList::find;
385
386 //- Find the edge within the face.
387 // \return the edge index or -1 if not found
388 label find(const Foam::edge& e) const;
389
390 //- Test the edge direction on the face
391 // \return
392 // - 0: edge not found on the face
393 // - +1: forward (counter-clockwise) on the face
394 // - -1: reverse (clockwise) on the face
395 int edgeDirection(const Foam::edge& e) const;
396
397 //- Find the longest edge on a face.
398 label longestEdge(const UList<point>& pts) const;
399
400
401 // Face splitting utilities
402
403 //- Number of triangles after splitting
404 inline label nTriangles() const noexcept;
405
406 //- Number of triangles after splitting
407 label nTriangles(const UList<point>& unused) const;
408
409 //- Split into triangles using existing points.
410 // Result in triFaces[triI..triI+nTri].
411 // Splits intelligently to maximize triangle quality.
412 // \Return number of faces created.
413 label triangles
414 (
415 const UList<point>& points,
416 label& triI,
417 faceList& triFaces
418 ) const;
419
420 //- Split into triangles using existing points.
421 // Append to DynamicList.
422 // \Return number of faces created.
423 template<int SizeMin>
424 label triangles
425 (
426 const UList<point>& points,
427 DynamicList<face, SizeMin>& triFaces
428 ) const;
429
430 //- Number of triangles and quads after splitting
431 // Returns the sum of both
432 label nTrianglesQuads
433 (
434 const UList<point>& points,
435 label& nTris,
436 label& nQuads
437 ) const;
438
439 //- Split into triangles and quads.
440 // Results in triFaces (starting at triI) and quadFaces
441 // (starting at quadI).
442 // Returns number of new faces created.
443 label trianglesQuads
444 (
445 const UList<point>& points,
446 label& triI,
447 label& quadI,
448 faceList& triFaces,
449 faceList& quadFaces
450 ) const;
451
452
453 //- True if the face has at least one vertex in common with other
454 inline bool connected(const labelUList& other) const;
455
456 //- True if the face has at least one vertex in common with other
457 template<unsigned N>
458 inline bool connected(const FixedList<label, N>& other) const;
459
460 //- Compare faces
461 // \return
462 // - 0: different
463 // - +1: identical
464 // - -1: same face, but different orientation
465 static int compare(const face& a, const face& b);
466
467 //- True if the faces have all the same vertices
468 static bool sameVertices(const face& a, const face& b);
469
470
471 // Member Operators
472
473 //- Increment (offset) vertices by given amount
474 inline void operator+=(const label vertexOffset);
475
476
477 // Hashing
478
479 //- The symmetric hash value for face.
480 // Starts with lowest vertex value and walks in the direction
481 // of the next lowest value.
482 static unsigned symmhash_code(const UList<label>& f, unsigned seed=0);
483
484 //- The hash value for face.
485 // Currently hashes as sequential contiguous data,
486 // but could/should be commutative
487 inline unsigned hash_code(unsigned seed=0) const
488 {
489 return Foam::Hasher(this->cdata(), this->size_bytes(), seed);
490 }
491
492 //- The symmetric hash value for face.
493 // Starts with lowest vertex value and walks in the direction
494 // of the next lowest value.
495 inline unsigned symmhash_code(unsigned seed=0) const
496 {
497 return face::symmhash_code(*this, seed);
498 }
499
500 //- Hashing functor for face
501 struct hasher
502 {
503 unsigned operator()(const face& obj, unsigned seed=0) const
504 {
505 return obj.hash_code(seed);
506 }
507 };
508
509 //- Symmetric hashing functor for face
510 struct symmHasher
511 {
512 unsigned operator()(const face& obj, unsigned seed=0) const
513 {
514 return face::symmhash_code(obj, seed);
515 }
516 };
517
518
519 // Housekeeping
520
521 //- Identical to edge()
522 Foam::edge faceEdge(label edgei) const { return this->edge(edgei); }
523};
524
525
526// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
527
528//- Hash face as contiguous (non-commutative) list data
529template<> struct Hash<face> : face::hasher {};
530
531
532//- Specialization to offset faces, used in ListListOps::combineOffset
533template<>
534struct offsetOp<face>
535{
536 face operator()(const face& x, const label offset) const
537 {
538 face f(x);
539 f += offset;
540 return f;
541 }
542};
543
544
545// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
546
547inline bool operator==(const face& a, const face& b);
548inline bool operator!=(const face& a, const face& b);
549
550
551// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
552
553} // End namespace Foam
554
555// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
556
557#include "faceI.H"
558
559#ifdef NoRepository
560 #include "faceTemplates.C"
561#endif
562
563// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
564
565#endif
566
567// ************************************************************************* //
label n
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-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:68
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
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:89
const label * cdata() const noexcept
std::streamsize size_bytes() const noexcept
static const Form zero
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
label trianglesQuads(const UList< point > &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split into triangles and quads.
Definition face.C:873
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
Return potential intersection with face with a ray starting at p, direction n (does not need to be no...
pointHit nearestPointClassify(const point &p, const UList< point > &meshPoints, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition faceI.H:158
void flip()
Flip the face in-place.
Definition face.C:492
static unsigned symmhash_code(const UList< label > &f, unsigned seed=0)
The symmetric hash value for face.
Definition face.C:414
pointHit intersection(const point &p, const vector &q, const point &ctr, const UList< point > &meshPoints, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Pair< point > box(const UList< point > &pts) const
The enclosing (bounding) box for the face.
Definition faceI.H:118
bool contains(const Foam::edge &e) const
True if face contains(edge).
Definition faceI.H:253
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:635
label nEdges() const noexcept
Return number of edges.
Definition faceI.H:151
label triangles(const UList< point > &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition face.C:845
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition faceI.H:105
unsigned hash_code(unsigned seed=0) const
The hash value for face.
Definition face.H:641
label nTrianglesQuads(const UList< point > &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition face.C:859
static int compare(const face &a, const face &b)
Compare faces.
Definition face.C:273
tensor inertia(const UList< point > &p, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference point and density specification.
Definition face.C:714
scalar areaInContact(const UList< point > &meshPoints, const scalarField &v) const
Return area in contact, given the displacement in vertices.
label find(const Foam::edge &e) const
Find the edge within the face.
Definition face.C:791
bool connected(const labelUList &other) const
True if the face has at least one vertex in common with other.
Definition faceI.H:226
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition face.C:816
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal(const UList< point > &p) const
Legacy name for areaNormal().
Definition face.H:251
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition faceI.H:174
label longestEdge(const UList< point > &pts) const
Find the longest edge on a face.
Definition face.C:885
proxType
The proximity classifications.
Definition face.H:128
@ POINT
Close to point.
Definition face.H:130
@ EDGE
Close to edge.
Definition face.H:131
@ NONE
Unknown proximity.
Definition face.H:129
pointHit nearestPoint(const point &p, const UList< point > &meshPoints) const
Return nearest point to face.
point average(const UList< point > &points) const
Average of face points: a quick estimate of the face centre.
Definition faceI.H:133
point centre(const UList< point > &points) const
Centre point of face.
Definition face.C:506
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition faceI.H:80
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
The sign for the side of the face plane the point is on, using three evenly distributed face points f...
label which(const label vertex) const
Find local vertex on face for the vertex label, same as find().
Definition faceI.H:196
scalar magSqr(const UList< point > &p) const
Magnitude squared of face area.
Definition faceI.H:111
label nextLabel(const label i) const
Next vertex on face.
Definition faceI.H:208
vector unitNormal(const UList< point > &p) const
The unit normal.
Definition faceI.H:97
label thisLabel(const label i) const
The vertex on face - identical to operator[], but with naming similar to nextLabel(),...
Definition faceI.H:202
edgeList rcEdges() const
Return list of edges in reverse walk order.
Definition face.C:770
label collapse()
Collapse face by removing duplicate vertex labels.
Definition face.C:467
face reverseFace() const
Return face with reverse direction.
Definition face.C:612
Foam::edge faceEdge(label edgei) const
Identical to edge().
Definition face.H:685
static bool sameVertices(const face &a, const face &b)
True if the faces have all the same vertices.
Definition face.C:374
edgeList edges() const
Return list of edges in forward walk order.
Definition face.C:749
vector areaNormal(const UList< point > &p) const
The area normal - with magnitude equal to area of face.
Definition face.C:569
label nTriangles() const noexcept
Number of triangles after splitting.
Definition faceI.H:220
unsigned symmhash_code(unsigned seed=0) const
The symmetric hash value for face.
Definition face.H:652
constexpr face() noexcept=default
Default construct.
label prevLabel(const label i) const
Previous vertex on face.
Definition faceI.H:214
static const char *const typeName
Definition face.H:137
A class for managing temporary objects.
Definition tmp.H:75
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
volScalarField & p
Forwards for various types of face lists.
const pointField & points
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
List< label > labelList
A List of labels.
Definition List.H:62
unsigned Hasher(const void *data, size_t len, unsigned seed=0)
Bob Jenkins's 96-bit mixer hashing function (lookup3).
Definition Hasher.C:575
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
labelList f(nPoints)
const pointField & pts
std::vector< Triangle > triangles
volScalarField & b
volScalarField & e
Hash function class. The default definition is for primitives. Non-primitives used to hash entries on...
Definition Hash.H:48
Hashing functor for face.
Definition face.H:661
unsigned operator()(const face &obj, unsigned seed=0) const
Definition face.H:662
Symmetric hashing functor for face.
Definition face.H:672
unsigned operator()(const face &obj, unsigned seed=0) const
Definition face.H:673
face operator()(const face &x, const label offset) const
Definition face.H:703
Offset operator for ListListOps::combineOffset().
T operator()(const T &x, const label offset) const
const Vector< label > N(dict.get< Vector< label > >("N"))