Loading...
Searching...
No Matches
triFace.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-2023 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::triFace
29
30Description
31 A triangular face using a FixedList of labels corresponding to mesh
32 vertices.
33
34See also
35 Foam::face, Foam::triangle
36
37SourceFiles
38 triFaceI.H
39 triFaceTemplates.C
40
41\*---------------------------------------------------------------------------*/
42
43#ifndef Foam_triFace_H
44#define Foam_triFace_H
45
46#include "FixedList.H"
47#include "edgeList.H"
48#include "pointHit.H"
49#include "intersection.H"
50#include "pointField.H"
51#include "triangle.H"
52#include "ListListOps.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56namespace Foam
57{
58
59// Forward Declarations
60class face;
61class triFace;
62
63inline bool operator==(const triFace& a, const triFace& b);
64inline bool operator!=(const triFace& a, const triFace& b);
66/*---------------------------------------------------------------------------*\
67 Class triFace Declaration
68\*---------------------------------------------------------------------------*/
69
70class triFace
71:
72 public FixedList<label, 3>
73{
74public:
75
76 // Generated Methods
77
78 //- The front() accessor (from FixedList) has no purpose
79 void front() = delete;
80
81 //- The back() accessor (from FixedList) has no purpose
82 void back() = delete;
83
84
85 // Constructors
86
87 //- Default construct, with invalid vertex labels (-1)
88 inline triFace();
89
90 //- Construct from three vertex labels
91 inline triFace(const label p0, const label p1, const label p2) noexcept;
92
93 //- Construct from an initializer list of three vertex labels
94 inline explicit triFace(std::initializer_list<label> list);
95
96 //- Copy construct from a list of three vertex labels.
97 inline explicit triFace(const labelUList& list);
98
99 //- Copy construct from a subset of vertex labels
100 inline triFace
101 (
102 const labelUList& list,
103 const FixedList<label, 3>& indices
104 );
105
106 //- Construct from Istream
107 inline triFace(Istream& is);
108
109
110 // Member Functions
111
112 // Access
113
114 //- The first vertex
115 label a() const noexcept { return get<0>(); }
116
117 //- The second vertex
118 label b() const noexcept { return get<1>(); }
119
120 //- The third vertex
121 label c() const noexcept { return get<2>(); }
122
123 //- The first vertex
124 label& a() noexcept { return get<0>(); }
125
126 //- The second vertex
127 label& b() noexcept { return get<1>(); }
129 //- The third vertex
130 label& c() noexcept { return get<2>(); }
131
132
133 // Queries
134
135 //- True if vertices are unique and non-negative.
136 inline bool good() const noexcept;
137
139 // Other
140
141 //- 'Collapse' face by marking duplicate vertex labels.
142 // Duplicates vertex labels are marked with '-1'
143 // (the lower vertex is retained).
144 // Return the collapsed size.
145 inline label collapse();
146
147 //- Flip the face in-place.
148 // The starting vertex of the original and flipped face are identical.
149 inline void flip();
150
151 //- Return the points corresponding to this face
152 inline pointField points(const UList<point>& pts) const;
154 //- Return triangle as a face
155 inline face triFaceFace() const;
156
157 //- Return the triangle
158 inline triPointRef tri(const UList<point>& points) const;
159
160 //- Return centre (centroid)
161 inline point centre(const UList<point>& points) const;
162
163 //- Calculate average value at centroid of face
164 template<class Type>
165 Type average(const UList<point>& unused, const Field<Type>& fld) const;
166
167 //- The area normal - with magnitude equal to area of face
168 inline vector areaNormal(const UList<point>& points) const;
169
170 //- The unit normal
171 inline vector unitNormal(const UList<point>& points) const;
172
173 //- Legacy name for areaNormal()
174 // \deprecated(2018-06) Deprecated for new use
175 FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
176 vector normal(const UList<point>& points) const
177 {
178 return areaNormal(points); // Legacy definition
179 }
180
181 //- Magnitude of face area
182 inline scalar mag(const UList<point>& points) const;
183
184 //- Magnitude squared of face area
185 inline scalar magSqr(const UList<point>& points) const;
186
187 //- The enclosing (bounding) box for the face
188 inline Pair<point> box(const UList<point>& points) const;
189
190 //- Number of triangles after splitting == 1
191 inline label nTriangles() const noexcept;
192
193 //- Return face with reverse direction
194 // The starting vertex of the original and reverse face are identical.
195 inline triFace reverseFace() const;
196
197 //- Find local vertex on face for the vertex label, same as find()
198 // \return position in face (0,1,2) or -1 if not found.
199 inline label which(const label vertex) const;
200
201 //- Next vertex on face
202 inline label nextLabel(const label i) const;
203
204 //- Previous vertex on face
205 inline label prevLabel(const label i) const;
206
207 //- The vertex on face - identical to operator[], but with naming
208 //- similar to nextLabel(), prevLabel()
209 inline label thisLabel(const label i) const;
210
211 //- Return swept-volume from old-points to new-points
212 inline scalar sweptVol
213 (
214 const UList<point>& opts,
215 const UList<point>& npts
216 ) const;
217
218 //- Return the inertia tensor, with optional reference
219 // point and density specification
220 inline tensor inertia
221 (
222 const UList<point>& points,
223 const point& refPt = vector::zero,
224 scalar density = 1.0
225 ) const;
226
227 //- Return point intersection with a ray starting at p, in direction q.
228 inline pointHit ray
229 (
230 const point& p,
231 const vector& q,
232 const UList<point>& points,
233 const intersection::algorithm = intersection::FULL_RAY,
234 const intersection::direction dir = intersection::VECTOR
235 ) const;
236
237 //- Fast intersection with a ray.
239 (
240 const point& p,
241 const vector& q,
242 const UList<point>& points,
243 const intersection::algorithm alg,
244 const scalar tol = 0.0
245 ) const;
246
248 (
249 const point& p,
250 const vector& q,
251 const point& ctr,
252 const UList<point>& points,
253 const intersection::algorithm alg,
254 const scalar tol = 0.0
255 ) const;
256
257 //- Return nearest point to face
259 (
260 const point& p,
261 const UList<point>& points
262 ) const;
263
264
265 //- Return nearest point to face and classify it:
266 // + near point (nearType=POINT, nearLabel=0, 1, 2)
267 // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
268 // Note: edges are counted from starting vertex so
269 // e.g. edge n is from f[n] to f[0], where the face has n + 1
270 // points
272 (
273 const point& p,
274 const UList<point>& points,
275 label& nearType,
276 label& nearLabel
277 ) const;
278
279 //- The sign for which side of the face plane the point is on.
280 // Uses the supplied tolerance for rounding around zero.
281 // \return
282 // - 0: on plane
283 // - +1: front-side
284 // - -1: back-side
285 inline int sign
286 (
287 const point& p,
288 const UList<point>& points,
289 const scalar tol = SMALL
290 ) const;
291
292 //- Return number of edges == 3
293 inline label nEdges() const noexcept;
294
295 //- Return i-th face edge (forward walk order).
296 // The edge 0 is from [0] to [1]
297 // and edge 1 is from [1] to [2]
298 inline Foam::edge edge(const label edgei) const;
299
300 //- Return vector of i-th face edge (forward walk order).
301 // The edge 0 is from [0] to [1]
302 // and edge 1 is from [1] to [2]
303 inline vector edge(const label edgei, const UList<point>& pts) const;
304
305 //- Return i-th face edge in reverse walk order.
306 // The rcEdge 0 is from [0] to [n-1]
307 // and rcEdge 1 is from [n-1] to [n-2]
308 inline Foam::edge rcEdge(const label edgei) const;
309
310 //- Return vector of i-th face edge in reverse walk order.
311 // The rcEdge 0 is from [0] to [n-1]
312 // and rcEdge 1 is from [n-1] to [n-2]
313 inline vector rcEdge(const label edgei, const UList<point>& pts) const;
314
315 //- Return list of edges in forward walk order.
316 // The edge 0 is from [0] to [1]
317 // and edge 1 is from [1] to [2]
318 inline edgeList edges() const;
319
320 //- Return list of edges in reverse walk order.
321 // The rcEdge 0 is from [0] to [n-1]
322 // and rcEdge 1 is from [n-1] to [n-2]
323 inline edgeList rcEdges() const;
324
325
326 // Searching
327
328 //- Regular contains(pointi) tests
329 using FixedList<label, 3>::contains;
330
331 //- True if face contains(edge)
332 inline bool contains(const Foam::edge& e) const;
333
334 //- Regular find pointi within face
335 using FixedList<label, 3>::find;
336
337 //- Find the edge within the face.
338 // \return the edge index or -1 if not found
339 inline label find(const Foam::edge& e) const;
340
341 //- Test the edge direction on the face
342 // \return
343 // - +1: forward (counter-clockwise) on the face
344 // - -1: reverse (clockwise) on the face
345 // - 0: edge not found on the face
346 inline int edgeDirection(const Foam::edge& e) const;
347
348 //- Compare triFaces
349 // \return:
350 // - 0: different
351 // - +1: identical
352 // - -1: same face, but different orientation
353 static inline int compare(const triFace& a, const triFace& b);
354
355
356 // Member Operators
357
358 //- Increment (offset) vertices by given amount
359 inline void operator+=(const label vertexOffset);
360
361
362 // Hashing
363
364 //- The (commutative) hash value for triFace
365 inline unsigned hash_code(unsigned seed=0) const
366 {
367 // Fortunately we don't need this very often
368 const uLabel t0((*this)[0]);
369 const uLabel t1((*this)[1]);
370 const uLabel t2((*this)[2]);
371
372 const uLabel val(t0*t1*t2 + t0+t1+t2);
373
374 return Foam::Hash<uLabel>()(val, seed);
375 }
376
377 //- Hashing functor for triFace (commutative)
378 struct hasher
379 {
380 unsigned operator()(const triFace& obj, unsigned seed=0) const
381 {
382 return obj.hash_code(seed);
383 }
384 };
385
386 //- Deprecated(2021-04) hashing functor. Use hasher()
387 // \deprecated(2021-04) - use hasher() functor
388 template<class Unused=bool>
389 struct Hash : triFace::hasher
390 {
391 FOAM_DEPRECATED_FOR(2021-04, "hasher()") Hash() {}
392 };
393
394
395 // Housekeeping
396
397 //- Same as good()
398 bool valid() const noexcept { return good(); }
399
400 //- Identical to edge()
401 Foam::edge faceEdge(label edgei) const { return this->edge(edgei); }
402};
403
404
405// * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
406
407//- Contiguous data for triFace
408template<> struct is_contiguous<triFace> : std::true_type {};
409
410//- Contiguous label data for triFace
411template<> struct is_contiguous_label<triFace> : std::true_type {};
412
413//- Hashing for triFace uses commutative (incremental) hash
414template<> struct Hash<triFace> : triFace::hasher {};
415
416
417//- Specialization to offset faces, used in ListListOps::combineOffset
418template<>
419struct offsetOp<triFace>
420{
421 triFace operator()(const triFace& x, const label offset) const
422 {
423 triFace f(x);
424 f += offset;
425 return f;
426 }
427};
428
429
430// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
431
432inline bool operator==(const triFace& a, const triFace& b);
433inline bool operator!=(const triFace& a, const triFace& b);
434
435
436// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
437
438} // End namespace Foam
439
440// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
441
442#include "triFaceI.H"
443
444#ifdef NoRepository
445 #include "triFaceTemplates.C"
446#endif
447
448// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449
450#endif
451
452// ************************************************************************* //
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))
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
label & get() noexcept
Definition FixedListI.H:152
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
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
Foam::intersection.
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
scalar sweptVol(const UList< point > &opts, const UList< point > &npts) const
Return swept-volume from old-points to new-points.
Definition triFaceI.H:248
vector unitNormal(const UList< point > &points) const
The unit normal.
Definition triFaceI.H:185
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal(const UList< point > &points) const
Legacy name for areaNormal().
Definition triFace.H:223
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition triFaceI.H:372
void flip()
Flip the face in-place.
Definition triFaceI.H:143
label b() const noexcept
The second vertex.
Definition triFace.H:133
pointHit ray(const point &p, const vector &q, const UList< point > &points, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return point intersection with a ray starting at p, in direction q.
Definition triFaceI.H:293
pointHit nearestPointClassify(const point &p, const UList< point > &points, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Definition triFaceI.H:344
void back()=delete
The back() accessor (from FixedList) has no purpose.
void front()=delete
The front() accessor (from FixedList) has no purpose.
pointHit nearestPoint(const point &p, const UList< point > &points) const
Return nearest point to face.
Definition triFaceI.H:334
scalar magSqr(const UList< point > &points) const
Magnitude squared of face area.
Definition triFaceI.H:197
pointHit intersection(const point &p, const vector &q, const UList< point > &points, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition triFaceI.H:307
bool contains(const Foam::edge &e) const
True if face contains(edge).
Definition triFaceI.H:488
label nEdges() const noexcept
Return number of edges == 3.
Definition triFaceI.H:366
unsigned hash_code(unsigned seed=0) const
The (commutative) hash value for triFace.
Definition triFace.H:489
triFace reverseFace() const
Return face with reverse direction.
Definition triFaceI.H:216
scalar mag(const UList< point > &points) const
Magnitude of face area.
Definition triFaceI.H:191
bool good() const noexcept
True if vertices are unique and non-negative.
Definition triFaceI.H:105
label find(const Foam::edge &e) const
Find the edge within the face.
Definition triFaceI.H:466
Pair< point > box(const UList< point > &points) const
The enclosing (bounding) box for the face.
Definition triFaceI.H:204
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition triFaceI.H:444
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition triFaceI.H:388
triFace()
Default construct, with invalid vertex labels (-1).
Definition triFaceI.H:56
label & b() noexcept
The second vertex.
Definition triFace.H:148
tensor inertia(const UList< point > &points, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition triFaceI.H:281
label c() const noexcept
The third vertex.
Definition triFace.H:138
point centre(const UList< point > &points) const
Return centre (centroid).
Definition triFaceI.H:173
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition triFaceI.H:149
int sign(const point &p, const UList< point > &points, const scalar tol=SMALL) const
The sign for which side of the face plane the point is on.
Definition triFaceI.H:356
static int compare(const triFace &a, const triFace &b)
Compare triFaces.
Definition triFaceI.H:27
label which(const label vertex) const
Find local vertex on face for the vertex label, same as find().
Definition triFaceI.H:223
label a() const noexcept
The first vertex.
Definition triFace.H:128
label nextLabel(const label i) const
Next vertex on face.
Definition triFaceI.H:235
bool valid() const noexcept
Same as good().
Definition triFace.H:529
label thisLabel(const label i) const
The vertex on face - identical to operator[], but with naming similar to nextLabel(),...
Definition triFaceI.H:229
edgeList rcEdges() const
Return list of edges in reverse walk order.
Definition triFaceI.H:427
label collapse()
'Collapse' face by marking duplicate vertex labels.
Definition triFaceI.H:116
face triFaceFace() const
Return triangle as a face.
Definition triFaceI.H:161
Foam::edge faceEdge(label edgei) const
Identical to edge().
Definition triFace.H:534
edgeList edges() const
Return list of edges in forward walk order.
Definition triFaceI.H:410
vector areaNormal(const UList< point > &points) const
The area normal - with magnitude equal to area of face.
Definition triFaceI.H:179
label nTriangles() const noexcept
Number of triangles after splitting == 1.
Definition triFaceI.H:210
Type average(const UList< point > &unused, const Field< Type > &fld) const
Calculate average value at centroid of face.
label & a() noexcept
The first vertex.
Definition triFace.H:143
label prevLabel(const label i) const
Previous vertex on face.
Definition triFaceI.H:241
label & c() noexcept
The third vertex.
Definition triFace.H:153
triPointRef tri(const UList< point > &points) const
Return the triangle.
Definition triFaceI.H:167
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
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
dimensionedScalar sign(const dimensionedScalar &ds)
Tensor< scalar > tensor
Definition symmTensor.H:57
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
face triFace(3)
labelList f(nPoints)
const pointField & pts
volScalarField & b
volScalarField & e
#define FOAM_DEPRECATED_FOR(since, replacement)
Definition stdFoam.H:43
Hash function class. The default definition is for primitives. Non-primitives used to hash entries on...
Definition Hash.H:48
A template class to specify if a data type is composed solely of Foam::label elements.
Definition contiguous.H:82
A template class to specify that a data type can be considered as being contiguous in memory.
Definition contiguous.H:70
triFace operator()(const triFace &x, const label offset) const
Definition triFace.H:562
Offset operator for ListListOps::combineOffset().
T operator()(const T &x, const label offset) const
Deprecated(2021-04) hashing functor. Use hasher().
Definition triFace.H:519
FOAM_DEPRECATED_FOR(2021-04, "hasher()") Hash()
Definition triFace.H:520
Hashing functor for triFace (commutative).
Definition triFace.H:505
unsigned operator()(const triFace &obj, unsigned seed=0) const
Definition triFace.H:506