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 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::triFace
29 
30 Description
31  A triangular face using a FixedList of labels corresponding to mesh
32  vertices.
33 
34 See also
35  Foam::face, Foam::triangle
36 
37 SourceFiles
38  triFaceI.H
39  triFaceTemplates.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef triFace_H
44 #define 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 
56 namespace Foam
57 {
58 
59 // Forward Declarations
60 class face;
61 class triFace;
62 
63 inline bool operator==(const triFace& a, const triFace& b);
64 inline bool operator!=(const triFace& a, const triFace& b);
65 
66 /*---------------------------------------------------------------------------*\
67  Class triFace Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class triFace
71 :
72  public FixedList<label, 3>
73 {
74 public:
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
112 
113  //- Return last (third) vertex label
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  // Other
124 
125  //- 'Collapse' face by marking duplicate point labels.
126  // Duplicates point labels are marked with '-1'
127  // (the lower vertex is retained).
128  // Return the collapsed size.
129  inline label collapse();
130 
131  //- Flip the face in-place.
132  // The starting points of the original and flipped face are identical.
133  inline void flip();
134 
135  //- Return the points corresponding to this face
136  inline pointField points(const UList<point>& pts) const;
137 
138  //- Return triangle as a face
139  inline face triFaceFace() const;
140 
141  //- Return the triangle
142  inline triPointRef tri(const UList<point>& points) const;
143 
144  //- Return centre (centroid)
145  inline point centre(const UList<point>& points) const;
146 
147  //- Calculate average value at centroid of face
148  template<class Type>
149  Type average(const UList<point>& unused, const Field<Type>& fld) const;
150 
151  //- The area normal - with magnitude equal to area of face
152  inline vector areaNormal(const UList<point>& points) const;
153 
154  //- The unit normal
155  inline vector unitNormal(const UList<point>& points) const;
156 
157  //- Legacy name for areaNormal()
158  // \deprecated(2018-06) Deprecated for new use
159  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
160  vector normal(const UList<point>& points) const
161  {
162  return areaNormal(points); // Legacy definition
163  }
164 
165  //- Magnitude of face area
166  inline scalar mag(const UList<point>& points) const;
167 
168  //- Number of triangles after splitting == 1
169  inline label nTriangles() const noexcept;
170 
171  //- Return face with reverse direction
172  // The starting points of the original and reverse face are identical.
173  inline triFace reverseFace() const;
174 
175  //- Find local index on face for the point label, same as find()
176  // \return position in face (0,1,2) or -1 if not found.
177  inline label which(const label pointLabel) const;
178 
179  //- Next vertex on face
180  inline label nextLabel(const label i) const;
181 
182  //- Previous vertex on face
183  inline label prevLabel(const label i) const;
184 
185  //- The vertex on face - identical to operator[], but with naming
186  //- similar to nextLabel(), prevLabel()
187  inline label thisLabel(const label i) const;
188 
189  //- Return swept-volume from old-points to new-points
190  inline scalar sweptVol
191  (
192  const UList<point>& opts,
193  const UList<point>& npts
194  ) const;
195 
196  //- Return the inertia tensor, with optional reference
197  // point and density specification
198  inline tensor inertia
199  (
200  const UList<point>& points,
201  const point& refPt = vector::zero,
202  scalar density = 1.0
203  ) const;
204 
205  //- Return point intersection with a ray starting at p, in direction q.
206  inline pointHit ray
207  (
208  const point& p,
209  const vector& q,
210  const UList<point>& points,
211  const intersection::algorithm = intersection::FULL_RAY,
212  const intersection::direction dir = intersection::VECTOR
213  ) const;
214 
215  //- Fast intersection with a ray.
216  inline pointHit intersection
217  (
218  const point& p,
219  const vector& q,
220  const UList<point>& points,
221  const intersection::algorithm alg,
222  const scalar tol = 0.0
223  ) const;
224 
225  inline pointHit intersection
226  (
227  const point& p,
228  const vector& q,
229  const point& ctr,
230  const UList<point>& points,
231  const intersection::algorithm alg,
232  const scalar tol = 0.0
233  ) const;
234 
235  //- Return nearest point to face
236  inline pointHit nearestPoint
237  (
238  const point& p,
239  const UList<point>& points
240  ) const;
241 
242 
243  //- Return nearest point to face and classify it:
244  // + near point (nearType=POINT, nearLabel=0, 1, 2)
245  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
246  // Note: edges are counted from starting vertex so
247  // e.g. edge n is from f[n] to f[0], where the face has n + 1
248  // points
250  (
251  const point& p,
252  const UList<point>& points,
253  label& nearType,
254  label& nearLabel
255  ) const;
256 
257  //- The sign for which side of the face plane the point is on.
258  // Uses the supplied tolerance for rounding around zero.
259  // \return
260  // - 0: on plane
261  // - +1: front-side
262  // - -1: back-side
263  inline int sign
264  (
265  const point& p,
266  const UList<point>& points,
267  const scalar tol = SMALL
268  ) const;
269 
270  //- Return number of edges == 3
271  inline label nEdges() const noexcept;
272 
273  //- Return i-th face edge (forward walk order).
274  // The edge 0 is from [0] to [1]
275  // and edge 1 is from [1] to [2]
276  inline Foam::edge edge(const label edgei) const;
277 
278  //- Return vector of i-th face edge (forward walk order).
279  // The edge 0 is from [0] to [1]
280  // and edge 1 is from [1] to [2]
281  inline vector edge(const label edgei, const UList<point>& pts) const;
282 
283  //- Return i-th face edge in reverse walk order.
284  // The rcEdge 0 is from [0] to [n-1]
285  // and rcEdge 1 is from [n-1] to [n-2]
286  inline Foam::edge rcEdge(const label edgei) const;
287 
288  //- Return vector of i-th face edge in reverse walk order.
289  // The rcEdge 0 is from [0] to [n-1]
290  // and rcEdge 1 is from [n-1] to [n-2]
291  inline vector rcEdge(const label edgei, const UList<point>& pts) const;
292 
293  //- Return list of edges in forward walk order.
294  // The edge 0 is from [0] to [1]
295  // and edge 1 is from [1] to [2]
296  inline edgeList edges() const;
297 
298  //- Return list of edges in reverse walk order.
299  // The rcEdge 0 is from [0] to [n-1]
300  // and rcEdge 1 is from [n-1] to [n-2]
301  inline edgeList rcEdges() const;
302 
303  //- Test the edge direction on the face
304  // \return
305  // - +1: forward (counter-clockwise) on the face
306  // - -1: reverse (clockwise) on the face
307  // - 0: edge not found on the face
308  inline int edgeDirection(const Foam::edge& e) const;
309 
310  //- Compare triFaces
311  // \return:
312  // - 0: different
313  // - +1: identical
314  // - -1: same face, but different orientation
315  static inline int compare(const triFace& a, const triFace& b);
316 
317 
318  // Hashing
319 
320  //- The (commutative) hash value for triFace
321  inline unsigned hash_code(unsigned seed=0) const
322  {
323  // Fortunately we don't need this very often
324  const uLabel t0((*this)[0]);
325  const uLabel t1((*this)[1]);
326  const uLabel t2((*this)[2]);
327 
328  const uLabel val(t0*t1*t2 + t0+t1+t2);
329 
330  return Foam::Hash<uLabel>()(val, seed);
331  }
332 
333  //- Hashing functor for triFace (commutative)
334  struct hasher
335  {
336  unsigned operator()(const triFace& obj, unsigned seed=0) const
337  {
338  return obj.hash_code(seed);
339  }
340  };
341 
342  //- Deprecated(2021-04) hashing functor. Use hasher()
343  // \deprecated(2021-04) - use hasher() functor
344  template<class Unused=bool>
345  struct Hash : triFace::hasher
346  {
347  FOAM_DEPRECATED_FOR(2021-04, "hasher()") Hash() {}
348  };
349 
350 
351  // Housekeeping
352 
353  //- Identical to edge()
354  Foam::edge faceEdge(label edgei) const { return this->edge(edgei); }
355 };
356 
357 
358 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
359 
360 //- Contiguous data for triFace
361 template<> struct is_contiguous<triFace> : std::true_type {};
362 
363 //- Contiguous label data for triFace
364 template<> struct is_contiguous_label<triFace> : std::true_type {};
365 
366 //- Hashing for triFace uses commutative (incremental) hash
367 template<> struct Hash<triFace> : triFace::hasher {};
368 
369 
370 //- Specialization to offset faces, used in ListListOps::combineOffset
371 template<>
372 struct offsetOp<triFace>
373 {
374  inline triFace operator()
375  (
376  const triFace& x,
377  const label offset
378  ) const
379  {
380  triFace result;
381 
382  forAll(x, i)
383  {
384  result[i] = x[i] + offset;
385  }
386  return result;
387  }
388 };
389 
390 
391 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
392 
393 inline bool operator==(const triFace& a, const triFace& b);
394 inline bool operator!=(const triFace& a, const triFace& b);
395 
396 
397 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
398 
399 } // End namespace Foam
400 
401 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
402 
403 #include "triFaceI.H"
404 
405 #ifdef NoRepository
406  #include "triFaceTemplates.C"
407 #endif
408 
409 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
410 
411 #endif
412 
413 // ************************************************************************* //
Foam::triFace::inertia
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:283
Foam::Tensor< scalar >
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::triFace::reverseFace
triFace reverseFace() const
Return face with reverse direction.
Definition: triFaceI.H:218
Foam::triFace::points
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition: triFaceI.H:147
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
intersection.H
triFaceI.H
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
Foam::triFace::edge
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition: triFaceI.H:374
Foam::triFace::thisLabel
label thisLabel(const label i) const
Definition: triFaceI.H:231
ListListOps.H
Foam::triFace::areaNormal
vector areaNormal(const UList< point > &points) const
The area normal - with magnitude equal to area of face.
Definition: triFaceI.H:187
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::triFace::triFaceFace
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:159
Foam::triFace::unitNormal
vector unitNormal(const UList< point > &points) const
The unit normal.
Definition: triFaceI.H:198
triPointRef.H
Foam::triFace::tri
triPointRef tri(const UList< point > &points) const
Return the triangle.
Definition: triFaceI.H:165
Foam::triFace::average
Type average(const UList< point > &unused, const Field< Type > &fld) const
Calculate average value at centroid of face.
Definition: triFaceTemplates.C:34
pointHit.H
Foam::triFace::flip
void flip()
Flip the face in-place.
Definition: triFaceI.H:141
Foam::triFace::collapse
label collapse()
'Collapse' face by marking duplicate point labels.
Definition: triFaceI.H:114
Foam::is_contiguous_label
A template class to specify if a data type is composed solely of Foam::label elements.
Definition: contiguous.H:83
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::triFace::rcEdges
edgeList rcEdges() const
Return list of edges in reverse walk order.
Definition: triFaceI.H:429
Foam::triFace::FOAM_DEPRECATED_FOR
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal(const UList< point > &points) const
Legacy name for areaNormal()
Definition: triFace.H:158
Foam::triFace::faceEdge
Foam::edge faceEdge(label edgei) const
Identical to edge()
Definition: triFace.H:353
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
Foam::triFace::second
label second() const
Return second vertex label.
Definition: triFace.H:119
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
Foam::triFace::which
label which(const label pointLabel) const
Find local index on face for the point label, same as find()
Definition: triFaceI.H:225
Foam::Hash
Hash function class. The default definition is for primitives. Non-primitives used to hash entries on...
Definition: Hash.H:53
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::triFace::nTriangles
label nTriangles() const noexcept
Number of triangles after splitting == 1.
Definition: triFaceI.H:212
Foam::triFace::sign
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:358
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::triFace::edges
edgeList edges() const
Return list of edges in forward walk order.
Definition: triFaceI.H:412
Foam::triFace::hash_code
unsigned hash_code(unsigned seed=0) const
The (commutative) hash value for triFace.
Definition: triFace.H:320
Foam::offsetOp
Offset operator for ListListOps::combineOffset()
Definition: ListListOps.H:101
Foam::triFace::Hash
Deprecated(2021-04) hashing functor. Use hasher()
Definition: triFace.H:344
Foam::triFace::sweptVol
scalar sweptVol(const UList< point > &opts, const UList< point > &npts) const
Return swept-volume from old-points to new-points.
Definition: triFaceI.H:250
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::triFace::centre
point centre(const UList< point > &points) const
Return centre (centroid)
Definition: triFaceI.H:176
edgeList.H
Foam::pointHit
PointHit< point > pointHit
A PointIndexHit for 3D points.
Definition: pointHit.H:44
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::FixedList< label, 3 >::operator[]
label & operator[](const label i)
Return element of FixedList.
Definition: FixedListI.H:412
Foam::triFace::edgeDirection
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition: triFaceI.H:446
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::triFace::nearestPoint
pointHit nearestPoint(const point &p, const UList< point > &points) const
Return nearest point to face.
Definition: triFaceI.H:336
Foam::triFace::nEdges
label nEdges() const noexcept
Return number of edges == 3.
Definition: triFaceI.H:368
pointField.H
Foam::triFace::hasher
Hashing functor for triFace (commutative)
Definition: triFace.H:333
Foam::triFace::nearestPointClassify
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:346
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
Foam::triFace::Hash::FOAM_DEPRECATED_FOR
FOAM_DEPRECATED_FOR(2021-04, "hasher()") Hash()
Definition: triFace.H:346
Foam::triFace::rcEdge
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition: triFaceI.H:390
Foam::triFace::ray
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:295
Foam::Vector< scalar >
Foam::triFace::mag
scalar mag(const UList< point > &points) const
Magnitude of face area.
Definition: triFaceI.H:206
Foam::List< edge >
Foam::triFace::second
label & second()
Return second vertex label.
Definition: triFace.H:116
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::triFace::intersection
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:309
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::intersection
Foam::intersection.
Definition: intersection.H:52
Foam::direction
uint8_t direction
Definition: direction.H:52
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::triFace::prevLabel
label prevLabel(const label i) const
Previous vertex on face.
Definition: triFaceI.H:243
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::triFace::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: triFaceI.H:237
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::triFace::hasher::operator()
unsigned operator()(const triFace &obj, unsigned seed=0) const
Definition: triFace.H:335
FixedList.H
triFaceTemplates.C
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::triFace::triFace
triFace()
Default construct, with invalid point labels (-1)
Definition: triFaceI.H:65
Foam::triFace::compare
static int compare(const triFace &a, const triFace &b)
Compare triFaces.
Definition: triFaceI.H:36
triFace
face triFace(3)
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:61
Foam::is_contiguous
A template class to specify that a data type can be considered as being contiguous in memory.
Definition: contiguous.H:75
Foam::zero
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:62