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-2020 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  //- 'Collapse' face by marking duplicate point labels.
109  // Duplicates point labels are marked with '-1'
110  // (the lower vertex is retained).
111  // Return the collapsed size.
112  inline label collapse();
113 
114  //- Flip the face in-place.
115  // The starting points of the original and flipped face are identical.
116  inline void flip();
117 
118  //- Return the points corresponding to this face
119  inline pointField points(const UList<point>& points) const;
120 
121  //- Return triangle as a face
122  inline face triFaceFace() const;
123 
124  //- Return the triangle
125  inline triPointRef tri(const UList<point>& points) const;
126 
127  //- Return centre (centroid)
128  inline point centre(const UList<point>& points) const;
129 
130  //- Calculate average value at centroid of face
131  template<class Type>
132  Type average(const UList<point>& unused, const Field<Type>& fld) const;
133 
134  //- The area normal - with magnitude equal to area of face
135  inline vector areaNormal(const UList<point>& points) const;
136 
137  //- The unit normal
138  inline vector unitNormal(const UList<point>& points) const;
139 
140  //- Legacy name for areaNormal()
141  // \deprecated(2018-06) Deprecated for new use
142  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
143  vector normal(const UList<point>& points) const
144  {
145  return areaNormal(points); // Legacy definition
146  }
147 
148  //- Magnitude of face area
149  inline scalar mag(const UList<point>& points) const;
150 
151  //- Number of triangles after splitting
152  inline label nTriangles() const;
153 
154  //- Return face with reverse direction
155  // The starting points of the original and reverse face are identical.
156  inline triFace reverseFace() const;
157 
158  //- Return true if the point label is found in face.
159  inline bool found(const label pointLabel) const;
160 
161  //- Find local index on face for the point label.
162  // \return position in face (0,1,2) or -1 if not found.
163  inline label which(const label pointLabel) const;
164 
165  //- Return swept-volume from old-points to new-points
166  inline scalar sweptVol
167  (
168  const UList<point>& opts,
169  const UList<point>& npts
170  ) const;
171 
172  //- Return the inertia tensor, with optional reference
173  // point and density specification
174  inline tensor inertia
175  (
176  const UList<point>& points,
177  const point& refPt = vector::zero,
178  scalar density = 1.0
179  ) const;
180 
181  //- Return point intersection with a ray starting at p, in direction q.
182  inline pointHit ray
183  (
184  const point& p,
185  const vector& q,
186  const UList<point>& points,
189  ) const;
190 
191  //- Fast intersection with a ray.
192  inline pointHit intersection
193  (
194  const point& p,
195  const vector& q,
196  const UList<point>& points,
197  const intersection::algorithm alg,
198  const scalar tol = 0.0
199  ) const;
200 
201  inline pointHit intersection
202  (
203  const point& p,
204  const vector& q,
205  const point& ctr,
206  const UList<point>& points,
207  const intersection::algorithm alg,
208  const scalar tol = 0.0
209  ) const;
210 
211  //- Return nearest point to face
212  inline pointHit nearestPoint
213  (
214  const point& p,
215  const UList<point>& points
216  ) const;
217 
218 
219  //- Return nearest point to face and classify it:
220  // + near point (nearType=POINT, nearLabel=0, 1, 2)
221  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
222  // Note: edges are counted from starting vertex so
223  // e.g. edge n is from f[n] to f[0], where the face has n + 1
224  // points
226  (
227  const point& p,
228  const UList<point>& points,
229  label& nearType,
230  label& nearLabel
231  ) const;
232 
233  //- The sign for which side of the face plane the point is on.
234  // Uses the supplied tolerance for rounding around zero.
235  // \return
236  // - 0: on plane
237  // - +1: front-side
238  // - -1: back-side
239  inline int sign
240  (
241  const point& p,
242  const UList<point>& points,
243  const scalar tol = SMALL
244  ) const;
245 
246  //- Return number of edges
247  inline label nEdges() const;
248 
249  //- Return edges in face point ordering,
250  // i.e. edges()[0] is edge between [0] and [1]
251  inline edgeList edges() const;
252 
253  //- Return n-th face edge
254  inline edge faceEdge(const label n) const;
255 
256  //- Return the edge direction on the face
257  // \return
258  // - +1: forward (counter-clockwise) on the face
259  // - -1: reverse (clockwise) on the face
260  // - 0: edge not found on the face
261  inline int edgeDirection(const edge& e) const;
262 
263  //- Compare triFaces
264  // \return:
265  // - 0: different
266  // - +1: identical
267  // - -1: same face, but different orientation
268  static inline int compare(const triFace& a, const triFace& b);
269 
270 
271  // Hashing
272 
273  //- The (commutative) hash-value for triFace
274  inline unsigned hashval(unsigned seed=0) const
275  {
276  // Fortunately we don't need this very often
277  const uLabel t0((*this)[0]);
278  const uLabel t1((*this)[1]);
279  const uLabel t2((*this)[2]);
280 
281  const uLabel val = (t0*t1*t2 + t0+t1+t2);
282 
283  return Foam::Hash<uLabel>()(val, seed);
284  }
285 
286  //- Hashing function class for triFace (commutative)
287  // Also useful for inheritance in sub-classes
288  template<class HashT=Foam::Hash<label>>
289  struct Hash
290  {
291  inline unsigned operator()
292  (
293  const triFace& obj,
294  unsigned seed=0
295  ) const
296  {
297  return obj.hashval(seed);
298  }
299  };
300 };
301 
302 
303 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
304 
305 //- Contiguous data for triFace
306 template<> struct is_contiguous<triFace> : std::true_type {};
307 
308 //- Contiguous label data for triFace
309 template<> struct is_contiguous_label<triFace> : std::true_type {};
310 
311 
312 //- Hash specialization for triFace as a commutative hash value.
313 template<>
314 struct Hash<triFace>
315 {
316  inline unsigned operator()(const triFace& obj, unsigned seed=0) const
317  {
318  return obj.hashval(seed);
319  }
320 };
321 
322 
323 //- Specialization to offset faces, used in ListListOps::combineOffset
324 template<>
325 struct offsetOp<triFace>
326 {
327  inline triFace operator()
328  (
329  const triFace& x,
330  const label offset
331  ) const
332  {
333  triFace result;
334 
335  forAll(x, i)
336  {
337  result[i] = x[i] + offset;
338  }
339  return result;
340  }
341 };
342 
343 
344 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
345 
346 inline bool operator==(const triFace& a, const triFace& b);
347 inline bool operator!=(const triFace& a, const triFace& b);
348 
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 } // End namespace Foam
353 
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355 
356 #include "triFaceI.H"
357 
358 #ifdef NoRepository
359  #include "triFaceTemplates.C"
360 #endif
361 
362 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
363 
364 #endif
365 
366 // ************************************************************************* //
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:271
Foam::intersection::direction
direction
Definition: intersection.H:66
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
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
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::intersection::FULL_RAY
Definition: intersection.H:74
Foam::triFace::found
bool found(const label pointLabel) const
Return true if the point label is found in face.
Definition: triFaceI.H:225
Foam::triFace::triFaceFace
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:159
Foam::Hash::operator()
unsigned operator()(const T &obj, unsigned seed=0) const
Definition: Hash.H:57
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::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:141
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:235
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:61
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triFace::which
label which(const label pointLabel) const
Find local index on face for the point label.
Definition: triFaceI.H:231
Foam::Hash
Hash function class. The default definition is for primitives, non-primitives used to hash entries on...
Definition: Hash.H:55
Foam::triFace::points
pointField points(const UList< point > &points) const
Return the points corresponding to this face.
Definition: triFaceI.H:147
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::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:346
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::Hash< triFace >::operator()
unsigned operator()(const triFace &obj, unsigned seed=0) const
Definition: triFace.H:315
Foam::triFace::edges
edgeList edges() const
Return edges in face point ordering,.
Definition: triFaceI.H:362
Foam::intersection::algorithm
algorithm
Definition: intersection.H:72
Foam::triFace::faceEdge
edge faceEdge(const label n) const
Return n-th face edge.
Definition: triFaceI.H:379
Foam::offsetOp
Offset operator for ListListOps::combineOffset()
Definition: ListListOps.H:101
Foam::triFace::Hash
Hashing function class for triFace (commutative)
Definition: triFace.H:288
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:238
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
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::intersection::VECTOR
Definition: intersection.H:68
Foam::triFace::nearestPoint
pointHit nearestPoint(const point &p, const UList< point > &points) const
Return nearest point to face.
Definition: triFaceI.H:324
pointField.H
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:334
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
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:283
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::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:297
Foam::UList< label >
Foam::triFace::hashval
unsigned hashval(unsigned seed=0) const
The (commutative) hash-value for triFace.
Definition: triFace.H:273
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::triFace::edgeDirection
int edgeDirection(const edge &e) const
Return the edge direction on the face.
Definition: triFaceI.H:389
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::zero
static const Vector< Cmpt > zero
Definition: VectorSpace.H:115
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
FixedList.H
triFaceTemplates.C
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
Foam::triFace::nTriangles
label nTriangles() const
Number of triangles after splitting.
Definition: triFaceI.H:212
Foam::triFace::nEdges
label nEdges() const
Return number of edges.
Definition: triFaceI.H:356
triFace
face triFace(3)
Foam::is_contiguous
A template class to specify that a data type can be considered as being contiguous in memory.
Definition: contiguous.H:75