triangle.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-2017 OpenFOAM Foundation
9  Copyright (C) 2018-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::triangle
29 
30 Description
31  A triangle primitive used to calculate face normals and swept volumes.
32 
33 SourceFiles
34  triangleI.H
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef triangle_H
39 #define triangle_H
40 
41 #include "intersection.H"
42 #include "vector.H"
43 #include "tensor.H"
44 #include "pointHit.H"
45 #include "Random.H"
46 #include "FixedList.H"
47 #include "UList.H"
48 #include "linePointRef.H"
49 #include "barycentric2D.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // Forward Declarations
57 class Istream;
58 class Ostream;
59 class plane;
60 class triPoints;
61 
62 template<class Point, class PointRef> class triangle;
63 
64 template<class Point, class PointRef>
66 
67 template<class Point, class PointRef>
69 
70 
71 // Common Typedefs
72 
73 //- A triangle using referred points
75 
76 
77 /*---------------------------------------------------------------------------*\
78  Class triangle Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 template<class Point, class PointRef>
82 class triangle
83 {
84 public:
85 
86  // Public Typedefs
87 
88  //- The point type
89  typedef Point point_type;
90 
91  //- Storage type for triangles originating from intersecting triangle
92  //- with another triangle
94 
95  //- The proximity classifications
96  enum proxType
97  {
98  NONE = 0,
100  EDGE
101  };
102 
103 
104  // Public classes
105 
106  // Classes for use in sliceWithPlane. What to do with decomposition
107  // of triangle.
108 
109  //- Dummy
110  class dummyOp
111  {
112  public:
113  inline void operator()(const triPoints&);
114  };
115 
116  //- Sum resulting areas
117  class sumAreaOp
118  {
119  public:
120  scalar area_;
121 
122  inline sumAreaOp();
123 
124  inline void operator()(const triPoints&);
125  };
126 
127  //- Store resulting tris
128  class storeOp
129  {
130  public:
132  label& nTris_;
133 
134  inline storeOp(triIntersectionList&, label&);
135 
136  inline void operator()(const triPoints&);
137  };
138 
139 
140 private:
141 
142  // Private Data
143 
144  PointRef a_, b_, c_;
145 
146 
147  // Private Member Functions
148 
149  //- Helper: calculate intersection point
150  inline static point planeIntersection
151  (
152  const FixedList<scalar, 3>& d,
153  const triPoints& t,
154  const label negI,
155  const label posI
156  );
157 
158  //- Helper: slice triangle with plane
159  template<class AboveOp, class BelowOp>
160  inline static void triSliceWithPlane
161  (
162  const plane& pln,
163  const triPoints& tri,
164  AboveOp& aboveOp,
165  BelowOp& belowOp
166  );
167 
168 
169 public:
170 
171  // Constructors
172 
173  //- Construct from three points
174  inline triangle(const Point& a, const Point& b, const Point& c);
175 
176  //- Construct from three points
177  inline triangle(const FixedList<Point, 3>& tri);
178 
179  //- Construct from three points in the list of points
180  // The indices could be from triFace etc.
181  inline triangle
182  (
183  const UList<Point>& points,
184  const FixedList<label, 3>& indices
185  );
186 
187  //- Construct from Istream
188  inline explicit triangle(Istream& is);
189 
190 
191  // Member Functions
192 
193  // Access
194 
195  //- Return first vertex
196  inline const Point& a() const;
197 
198  //- Return second vertex
199  inline const Point& b() const;
200 
201  //- Return third vertex
202  inline const Point& c() const;
203 
204 
205  // Properties
206 
207  //- Return centre (centroid)
208  inline Point centre() const;
209 
210  //- The area normal - with magnitude equal to area of triangle
211  inline vector areaNormal() const;
212 
213  //- Return unit normal
214  inline vector unitNormal() const;
215 
216  //- Legacy name for areaNormal().
217  // \deprecated(2018-06) Deprecated for new use
218  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
219  vector normal() const
220  {
221  return areaNormal();
222  }
223 
224  //- Return scalar magnitude
225  inline scalar mag() const;
226 
227  //- Return circum-centre
228  inline Point circumCentre() const;
229 
230  //- Return circum-radius
231  inline scalar circumRadius() const;
232 
233  //- Return quality: Ratio of triangle and circum-circle
234  // area, scaled so that an equilateral triangle has a
235  // quality of 1
236  inline scalar quality() const;
237 
238  //- Return swept-volume
239  inline scalar sweptVol(const triangle& t) const;
240 
241  //- Return the inertia tensor, with optional reference
242  // point and density specification
243  inline tensor inertia
244  (
245  PointRef refPt = Zero,
246  scalar density = 1.0
247  ) const;
248 
249  //- Return a random point on the triangle from a uniform
250  // distribution
251  inline Point randomPoint(Random& rndGen) const;
252 
253  //- Calculate the point from the given barycentric coordinates.
254  inline Point barycentricToPoint(const barycentric2D& bary) const;
255 
256  //- Calculate the barycentric coordinates from the given point
257  inline barycentric2D pointToBarycentric(const point& pt) const;
258 
259  //- Calculate the barycentric coordinates from the given point.
260  // Returns the determinant.
261  inline scalar pointToBarycentric
262  (
263  const point& pt,
264  barycentric2D& bary
265  ) const;
266 
267  //- Return point intersection with a ray.
268  // For a hit, the distance is signed. Positive number
269  // represents the point in front of triangle.
270  // In case of miss pointHit is set to nearest point
271  // on triangle and its distance to the distance between
272  // the original point and the plane intersection point
273  inline pointHit ray
274  (
275  const point& p,
276  const vector& q,
279  ) const;
280 
281  //- Fast intersection with a ray.
282  // For a hit, the pointHit.distance() is the line parameter t :
283  // intersection=p+t*q. Only defined for VISIBLE, FULL_RAY or
284  // HALF_RAY. tol increases the virtual size of the triangle
285  // by a relative factor.
286  inline pointHit intersection
287  (
288  const point& p,
289  const vector& q,
290  const intersection::algorithm alg,
291  const scalar tol = 0.0
292  ) const;
293 
294  //- Find the nearest point to p on the triangle and classify it:
295  // + near point (nearType=POINT, nearLabel=0, 1, 2)
296  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
297  // Note: edges are counted from starting
298  // vertex so e.g. edge 2 is from f[2] to f[0]
300  (
301  const point& p,
302  label& nearType,
303  label& nearLabel
304  ) const;
305 
306  //- Return nearest point to p on triangle
307  inline pointHit nearestPoint(const point& p) const;
308 
309  //- Classify nearest point to p in triangle plane
310  // w.r.t. triangle edges and points. Returns inside
311  // (true)/outside (false).
312  bool classify
313  (
314  const point& p,
315  label& nearType,
316  label& nearLabel
317  ) const;
318 
319  //- Return nearest point to line on triangle. Returns hit if
320  // point is inside triangle. Sets edgePoint to point on edge
321  // (hit if nearest is inside line)
322  inline pointHit nearestPoint
323  (
324  const linePointRef& edge,
325  pointHit& edgePoint
326  ) const;
327 
328  //- The sign for which side of the face plane the point is on.
329  // Uses the supplied tolerance for rounding around zero.
330  // \return
331  // - 0: on plane
332  // - +1: front-side
333  // - -1: back-side
334  inline int sign(const point& p, const scalar tol = SMALL) const;
335 
336  //- Decompose triangle into triangles above and below plane
337  template<class AboveOp, class BelowOp>
338  inline void sliceWithPlane
339  (
340  const plane& pln,
341  AboveOp& aboveOp,
342  BelowOp& belowOp
343  ) const;
344 
345  //- Decompose triangle into triangles inside and outside
346  // (with respect to user provided normal) other
347  // triangle.
348  template<class InsideOp, class OutsideOp>
349  inline void triangleOverlap
350  (
351  const vector& n,
352  const triangle<Point, PointRef>& tri,
353  InsideOp& insideOp,
354  OutsideOp& outsideOp
355  ) const;
356 
357 
358  // IOstream operators
359 
360  friend Istream& operator>> <Point, PointRef>
361  (
362  Istream&,
363  triangle&
364  );
365 
366  friend Ostream& operator<< <Point, PointRef>
367  (
368  Ostream&,
369  const triangle&
370  );
371 };
372 
373 
374 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375 
376 } // End namespace Foam
377 
378 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
379 
380 #include "triangleI.H"
381 
382 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
383 
384 #ifdef NoRepository
385 # include "triangle.C"
386 #endif
387 
388 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
389 
390 #endif
391 
392 // ************************************************************************* //
Foam::intersection::direction
direction
Definition: intersection.H:66
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::Tensor< scalar >
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::triangle::triIntersectionList
FixedList< triPoints, 27 > triIntersectionList
Definition: triangle.H:92
Foam::triangle::storeOp::operator()
void operator()(const triPoints &)
Definition: triangleI.H:873
intersection.H
Foam::triangle::storeOp
Store resulting tris.
Definition: triangle.H:127
Foam::triangle::randomPoint
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform.
Definition: triangleI.H:254
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::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::triangle::sumAreaOp::sumAreaOp
sumAreaOp()
Definition: triangleI.H:843
Foam::triangle::sign
int sign(const point &p, const scalar tol=SMALL) const
The sign for which side of the face plane the point is on.
Definition: triangleI.H:823
Foam::triangle::classify
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:685
Foam::triangle::inertia
tensor inertia(PointRef refPt=Zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:218
triangle.C
Foam::triangle::dummyOp
Dummy.
Definition: triangle.H:109
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::Barycentric2D
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant.
Definition: Barycentric2D.H:54
Foam::triangle::sumAreaOp::operator()
void operator()(const triPoints &)
Definition: triangleI.H:851
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
pointHit.H
tensor.H
Foam::triangle::centre
Point centre() const
Return centre (centroid)
Definition: triangleI.H:105
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::triangle::unitNormal
vector unitNormal() const
Return unit normal.
Definition: triangleI.H:119
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:61
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triPoints
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:52
Foam::triangle::c
const Point & c() const
Return third vertex.
Definition: triangleI.H:98
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::triangle::point_type
Point point_type
The point type.
Definition: triangle.H:88
Foam::triangle::a
const Point & a() const
Return first vertex.
Definition: triangleI.H:86
Foam::triangle::NONE
Unknown proximity.
Definition: triangle.H:97
Foam::intersection::algorithm
algorithm
Definition: intersection.H:72
Foam::triangle::ray
pointHit ray(const point &p, const vector &q, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return point intersection with a ray.
Definition: triangleI.H:323
Foam::triangle::sumAreaOp::area_
scalar area_
Definition: triangle.H:119
Foam::triangle::b
const Point & b() const
Return second vertex.
Definition: triangleI.H:92
Foam::triangle::circumRadius
scalar circumRadius() const
Return circum-radius.
Definition: triangleI.H:162
Foam::triangle::POINT
Close to point.
Definition: triangle.H:98
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::intersection::VECTOR
Definition: intersection.H:68
triangleI.H
Random.H
Foam::triangle::FOAM_DEPRECATED_FOR
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal() const
Legacy name for areaNormal().
Definition: triangle.H:217
Foam::triangle::sumAreaOp
Sum resulting areas.
Definition: triangle.H:116
Foam::triangle::mag
scalar mag() const
Return scalar magnitude.
Definition: triangleI.H:128
Foam::triangle::nearestPointClassify
pointHit nearestPointClassify(const point &p, label &nearType, label &nearLabel) const
Find the nearest point to p on the triangle and classify it:
Definition: triangleI.H:522
Foam::triangle::EDGE
Close to edge.
Definition: triangle.H:99
UList.H
barycentric2D.H
Foam::triangle::proxType
proxType
The proximity classifications.
Definition: triangle.H:95
Foam::Vector< scalar >
Foam::triangle::storeOp::nTris_
label & nTris_
Definition: triangle.H:131
Foam::triangle::quality
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition: triangleI.H:183
Foam::triangle::circumCentre
Point circumCentre() const
Return circum-centre.
Definition: triangleI.H:135
Foam::triangle::dummyOp::operator()
void operator()(const triPoints &)
Definition: triangleI.H:836
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::triangle::intersection
pointHit intersection(const point &p, const vector &q, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition: triangleI.H:439
Foam::line
A line primitive.
Definition: line.H:59
vector.H
Foam::triangle::triangle
triangle(const Point &a, const Point &b, const Point &c)
Construct from three points.
Definition: triangleI.H:38
Foam::triangle::storeOp::tris_
triIntersectionList & tris_
Definition: triangle.H:130
Point
CGAL::Point_3< K > Point
Definition: CGALIndexedPolyhedron.H:53
rndGen
Random rndGen
Definition: createFields.H:23
Foam::triangle::nearestPoint
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:671
Foam::triangle::pointToBarycentric
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: triangleI.H:272
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::triangle::barycentricToPoint
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition: triangleI.H:262
FixedList.H
Foam::triangle::areaNormal
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition: triangleI.H:112
Foam::triangle::sweptVol
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition: triangleI.H:199
Foam::triangle::sliceWithPlane
void sliceWithPlane(const plane &pln, AboveOp &aboveOp, BelowOp &belowOp) const
Decompose triangle into triangles above and below plane.
Definition: triangle.C:115
Foam::triPointRef
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:73
linePointRef.H
Foam::triangle::triangleOverlap
void triangleOverlap(const vector &n, const triangle< Point, PointRef > &tri, InsideOp &insideOp, OutsideOp &outsideOp) const
Decompose triangle into triangles inside and outside.
Definition: triangle.C:128
Foam::triangle::storeOp::storeOp
storeOp(triIntersectionList &, label &)
Definition: triangleI.H:861