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