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