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