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-2022 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  Uses referred points.
33 
34 SourceFiles
35  triangleI.H
36  triangle.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef Foam_triangle_H
41 #define Foam_triangle_H
42 
43 #include "triangleFwd.H"
44 #include "intersection.H"
45 #include "vector.H"
46 #include "tensor.H"
47 #include "pointHit.H"
48 #include "Random.H"
49 #include "FixedList.H"
50 #include "UList.H"
51 #include "line.H"
52 #include "Pair.H"
53 #include "Tuple2.H"
54 #include "barycentric2D.H"
55 #include "treeBoundBox.H"
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 namespace Foam
60 {
61 
62 // Forward Declarations
63 class plane;
64 
65 template<class Point, class PointRef>
66 inline Istream& operator>>(Istream&, triangle<Point, PointRef>&);
67 
68 template<class Point, class PointRef>
69 inline Ostream& operator<<(Ostream&, const triangle<Point, PointRef>&);
70 
71 
72 /*---------------------------------------------------------------------------*\
73  Class triPoints Declaration
74 \*---------------------------------------------------------------------------*/
75 
76 //- Triangle point storage. Default constructable (triangle is not)
77 class triPoints
78 :
79  public FixedList<point, 3>
80 {
81 public:
82 
83  // Generated Methods
84 
85  //- Default construct
86  triPoints() = default;
87 
88  //- The front() accessor (from FixedList) has no purpose
89  void front() = delete;
90 
91  //- The back() accessor (from FixedList) has no purpose
92  void back() = delete;
93 
94 
95  // Constructors
96 
97  //- Construct from three points
98  inline triPoints(const point& p0, const point& p1, const point& p2);
99 
100  //- Construct from point references
101  inline explicit triPoints(const triPointRef& pts);
102 
103  //- Construct from three points
104  inline triPoints(const FixedList<point, 3>& pts);
105 
106  //- Copy construct from subset of points
107  inline triPoints
108  (
109  const UList<point>& points,
110  const FixedList<label, 3>& indices
111  );
112 
113 
114  // Member Functions
115 
116  //- The first vertex
117  const point& a() const noexcept { return get<0>(); }
118 
119  //- The second vertex
120  const point& b() const noexcept { return get<1>(); }
121 
122  //- The third vertex
123  const point& c() const noexcept { return get<2>(); }
124 
125  //- The first vertex
126  point& a() noexcept { return get<0>(); }
127 
128  //- The second vertex
129  point& b() noexcept { return get<1>(); }
131  //- The third vertex
132  point& c() noexcept { return get<2>(); }
133 
134  //- Flip triangle orientation by swapping second and third vertices
135  inline void flip();
136 
137  //- Return as triangle reference
138  inline triPointRef tri() const;
139 
141  // Properties
142 
143  //- Return centre (centroid)
144  inline point centre() const;
146  //- The area normal - with magnitude equal to area of triangle
147  inline vector areaNormal() const;
148 
149  //- Return unit normal
150  inline vector unitNormal() const;
151 
152  //- The magnitude of the triangle area
153  inline scalar mag() const;
154 
155  //- The enclosing (bounding) box for the triangle
156  inline Pair<point> box() const;
157 
158  //- Edge vector opposite point a(): from b() to c()
159  inline vector vecA() const;
160 
161  //- Edge vector opposite point b(): from c() to a()
162  inline vector vecB() const;
163 
164  //- Edge vector opposite point c(): from a() to b()
165  inline vector vecC() const;
166 };
167 
168 
169 /*---------------------------------------------------------------------------*\
170  Class triangle Declaration
171 \*---------------------------------------------------------------------------*/
172 
173 template<class Point, class PointRef>
174 class triangle
175 {
176 public:
177 
178  // Public Typedefs
179 
180  //- The point type
181  typedef Point point_type;
182 
183  //- Storage type for triangles originating from intersecting triangle
184  //- with another triangle
186 
187  //- Proximity classifications
188  enum proxType
189  {
190  NONE = 0,
191  POINT,
192  EDGE
193  };
194 
195 
196  // Public Classes
197 
198  // Classes for use in sliceWithPlane.
199  // What to do with decomposition of triangle.
200 
201  //- Dummy
202  class dummyOp
203  {
204  public:
205  inline void operator()(const triPoints&);
206  };
207 
208  //- Sum resulting areas
209  class sumAreaOp
210  {
211  public:
212  scalar area_;
213 
214  inline sumAreaOp();
215 
216  inline void operator()(const triPoints&);
217  };
218 
219  //- Store resulting tris
220  class storeOp
221  {
222  public:
224  label& nTris_;
225 
226  inline storeOp(triIntersectionList&, label&);
227 
228  inline void operator()(const triPoints&);
229  };
230 
231 
232 private:
233 
234  // Private Data
235 
236  //- Reference to the first triangle point
237  PointRef a_;
238 
239  //- Reference to the second triangle point
240  PointRef b_;
242  //- Reference to the third triangle point
243  PointRef c_;
244 
245 
246  // Private Member Functions
247 
248  //- Helper: calculate intersection point
249  inline static point planeIntersection
250  (
251  const FixedList<scalar, 3>& d,
252  const triPoints& t,
253  const label negI,
254  const label posI
255  );
256 
257  //- Helper: slice triangle with plane
258  template<class AboveOp, class BelowOp>
259  inline static void triSliceWithPlane
260  (
261  const plane& pln,
262  const triPoints& tri,
263  AboveOp& aboveOp,
264  BelowOp& belowOp
265  );
266 
267 
268 public:
269 
270  // Constructors
271 
272  //- Construct from three points
273  inline triangle(const Point& p0, const Point& p1, const Point& p2);
274 
275  //- Construct from three points
276  inline triangle(const FixedList<Point, 3>& pts);
277 
278  //- Construct from three points in the list of points
279  // The indices could be from triFace etc.
280  inline triangle
281  (
282  const UList<Point>& points,
283  const FixedList<label, 3>& indices
284  );
285 
286  //- Construct from Istream
287  inline explicit triangle(Istream& is);
288 
289 
290  // Member Functions
291 
292  // Access
293 
294  //- The first vertex
295  const Point& a() const noexcept { return a_; }
296 
297  //- The second vertex
298  const Point& b() const noexcept { return b_; }
299 
300  //- The third vertex
301  const Point& c() const noexcept { return c_; }
302 
303 
304  // Triangle properties (static calculations)
305 
306  //- The centre (centroid) of three points
307  inline static Point centre
308  (
309  const Point& p0,
310  const Point& p1,
311  const Point& p2
312  );
313 
314  //- The area normal for a triangle defined by three points
315  //- (right-hand rule). Magnitude equal to area of the triangle
316  inline static vector areaNormal
317  (
318  const Point& p0,
319  const Point& p1,
320  const Point& p2
321  );
322 
323  //- The unit normal for a triangle defined by three points
324  //- (right-hand rule).
325  inline static vector unitNormal
326  (
327  const Point& p0,
328  const Point& p1,
329  const Point& p2
330  );
331 
332  //- The enclosing (bounding) box for three points
333  inline static Pair<Point> box
334  (
335  const Point& p0,
336  const Point& p1,
337  const Point& p2
338  );
339 
340 
341  // Properties
342 
343  //- Return centre (centroid)
344  inline Point centre() const;
345 
346  //- The area normal - with magnitude equal to area of triangle
347  inline vector areaNormal() const;
348 
349  //- Return unit normal
350  inline vector unitNormal() const;
351 
352  //- Legacy name for areaNormal().
353  // \deprecated(2018-06) Deprecated for new use
354  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
355  vector normal() const
356  {
357  return areaNormal();
358  }
359 
360  //- The magnitude of the triangle area
361  inline scalar mag() const;
362 
363  //- The enclosing (bounding) box for the triangle
364  inline Pair<Point> box() const;
365 
366  //- Edge vector opposite point a(): from b() to c()
367  inline Point vecA() const;
368 
369  //- Edge vector opposite point b(): from c() to a()
370  inline Point vecB() const;
372  //- Edge vector opposite point c(): from a() to b()
373  inline Point vecC() const;
374 
375 
376  // Properties
377 
378  //- Return circum-centre
379  inline Point circumCentre() const;
380 
381  //- Return circum-radius
382  inline scalar circumRadius() const;
383 
384  //- Return quality: Ratio of triangle and circum-circle
385  // area, scaled so that an equilateral triangle has a
386  // quality of 1
387  inline scalar quality() const;
388 
389  //- Return swept-volume
390  inline scalar sweptVol(const triangle& t) const;
391 
392  //- Return the inertia tensor, with optional reference
393  // point and density specification
394  inline tensor inertia
395  (
396  PointRef refPt = Zero,
397  scalar density = 1.0
398  ) const;
399 
400  //- Return a random point on the triangle from a uniform
401  //- distribution
402  inline Point randomPoint(Random& rndGen) const;
403 
404  //- Calculate the point from the given barycentric coordinates.
405  inline Point barycentricToPoint(const barycentric2D& bary) const;
406 
407  //- Calculate the barycentric coordinates from the given point
408  inline barycentric2D pointToBarycentric(const point& pt) const;
409 
410  //- Calculate the barycentric coordinates from the given point.
411  // Returns the determinant.
412  inline scalar pointToBarycentric
413  (
414  const point& pt,
415  barycentric2D& bary
416  ) const;
417 
418  //- Fast intersection detection with a plane.
419  inline bool intersects
420  (
421  const point& origin,
422  const vector& normal
423  ) const;
424 
425  //- Fast intersection detection with an \b axis plane.
426  inline bool intersects
427  (
429  const point& origin,
431  const vector::components axis
432  ) const;
433 
434  //- Return point intersection with a ray.
435  // For a hit, the distance is signed. Positive number
436  // represents the point in front of triangle.
437  // In case of miss pointHit is set to nearest point
438  // on triangle and its distance to the distance between
439  // the original point and the plane intersection point
440  inline pointHit ray
441  (
442  const point& p,
443  const vector& q,
446  ) const;
447 
448  //- Fast intersection with a ray.
449  // For a hit, the pointHit.distance() is the line parameter t :
450  // intersection=p+t*q. Only defined for VISIBLE, FULL_RAY or
451  // HALF_RAY. tol increases the virtual size of the triangle
452  // by a relative factor.
453  inline pointHit intersection
454  (
455  const point& p,
456  const vector& q,
457  const intersection::algorithm alg,
458  const scalar tol = 0.0
459  ) const;
460 
461  //- Find the nearest point to p on the triangle and classify it:
462  // + near point (nearType=POINT, nearLabel=0, 1, 2)
463  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
464  // Note: edges are counted from starting
465  // vertex so e.g. edge 2 is from f[2] to f[0]
467  (
468  const point& p,
469  label& nearType,
470  label& nearLabel
471  ) const;
472 
473  //- Return nearest point to p on triangle
474  inline pointHit nearestPoint(const point& p) const;
475 
476  //- Classify nearest point to p in triangle plane
477  // w.r.t. triangle edges and points. Returns inside
478  // (true)/outside (false).
479  bool classify
480  (
481  const point& p,
482  label& nearType,
483  label& nearLabel
484  ) const;
485 
486  //- Return nearest point to line on triangle. Returns hit if
487  // point is inside triangle. Sets edgePoint to point on edge
488  // (hit if nearest is inside line)
489  inline pointHit nearestPoint
490  (
491  const linePointRef& edge,
492  pointHit& edgePoint
493  ) const;
494 
495  //- The sign for which side of the face plane the point is on.
496  // Uses the supplied tolerance for rounding around zero.
497  // \return
498  // - 0: on plane
499  // - +1: above plane
500  // - -1: below plane
501  inline int sign(const point& p, const scalar tol = SMALL) const;
502 
503  //- Decompose triangle into triangles above and below plane
504  template<class AboveOp, class BelowOp>
505  inline void sliceWithPlane
506  (
507  const plane& pln,
508  AboveOp& aboveOp,
509  BelowOp& belowOp
510  ) const;
511 
512  //- Decompose triangle into triangles inside and outside
513  // (with respect to user provided normal) other
514  // triangle.
515  template<class InsideOp, class OutsideOp>
516  inline void triangleOverlap
517  (
518  const vector& n,
519  const triangle<Point, PointRef>& tri,
520  InsideOp& insideOp,
521  OutsideOp& outsideOp
522  ) const;
523 
524 
525  // IOstream Operators
526 
527  friend Istream& operator>> <Point, PointRef>
528  (
529  Istream&,
530  triangle&
531  );
532 
533  friend Ostream& operator<< <Point, PointRef>
534  (
535  Ostream&,
536  const triangle&
537  );
538 };
539 
540 
541 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
542 
543 } // End namespace Foam
544 
545 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
546 
547 #include "triangleI.H"
548 
549 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
550 
551 #ifdef NoRepository
552 # include "triangle.C"
553 #endif
554 
555 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
556 
557 #endif
558 
559 // ************************************************************************* //
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition: triangleI.H:351
const point & c() const noexcept
The third vertex.
Definition: triangle.H:140
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
A line primitive.
Definition: line.H:52
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: triangleI.H:424
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
scalar mag() const
The magnitude of the triangle area.
Definition: triangleI.H:280
Unknown proximity.
Definition: triangle.H:239
vector vecC() const
Edge vector opposite point c(): from a() to b()
Definition: triangleI.H:250
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition: triangleI.H:335
void sliceWithPlane(const plane &pln, AboveOp &aboveOp, BelowOp &belowOp) const
Decompose triangle into triangles above and below plane.
Definition: triangle.C:107
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Pair< Point > box() const
The enclosing (bounding) box for the triangle.
Definition: triangleI.H:207
vector vecB() const
Edge vector opposite point b(): from c() to a()
Definition: triangleI.H:244
tensor inertia(PointRef refPt=Zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:370
Random rndGen
Definition: createFields.H:23
const Point & a() const noexcept
The first vertex.
Definition: triangle.H:371
Point point_type
The point type.
Definition: triangle.H:226
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:712
triPointRef tri() const
Return as triangle reference.
Definition: triangleI.H:258
void flip()
Flip triangle orientation by swapping second and third vertices.
Definition: triangleI.H:264
triPoints()=default
Default construct.
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:90
const Point & c() const noexcept
The third vertex.
Definition: triangle.H:381
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:1011
Point vecB() const
Edge vector opposite point b(): from c() to a()
Definition: triangleI.H:226
Close to point.
Definition: triangle.H:240
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant...
Definition: Barycentric2D.H:50
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal() const
Legacy name for areaNormal().
Definition: triangle.H:451
void back()=delete
The back() accessor (from FixedList) has no purpose.
triangle(const Point &p0, const Point &p1, const Point &p2)
Construct from three points.
Definition: triangleI.H:67
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:875
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:861
proxType
Proximity classifications.
Definition: triangle.H:237
vector unitNormal() const
Return unit normal.
Definition: triangleI.H:178
Close to edge.
Definition: triangle.H:241
const pointField & points
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:515
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
const Point & b() const noexcept
The second vertex.
Definition: triangle.H:376
Point vecA() const
Edge vector opposite point a(): from b() to c()
Definition: triangleI.H:220
Istream & operator>>(Istream &, directionInfo &)
Triangle point storage. Default constructable (triangle is not)
Definition: triangle.H:74
FixedList< triPoints, 27 > triIntersectionList
Storage type for triangles originating from intersecting triangle with another triangle.
Definition: triangle.H:232
void operator()(const triPoints &)
Definition: triangleI.H:1061
Point centre() const
Return centre (centroid)
Definition: triangleI.H:126
void operator()(const triPoints &)
Definition: triangleI.H:1024
scalar mag() const
The magnitude of the triangle area.
Definition: triangleI.H:273
Vector< scalar > vector
Definition: vector.H:57
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
scalar circumRadius() const
Return circum-radius.
Definition: triangleI.H:314
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:105
Random number generator.
Definition: Random.H:55
storeOp(triIntersectionList &, label &)
Definition: triangleI.H:1049
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
void front()=delete
The front() accessor (from FixedList) has no purpose.
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition: triangleI.H:157
const point & a() const noexcept
The first vertex.
Definition: triangle.H:130
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition: triangleI.H:414
bool intersects(const point &origin, const vector &normal) const
Fast intersection detection with a plane.
Definition: triangleI.H:475
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:120
CGAL::Point_3< K > Point
Pair< point > box() const
The enclosing (bounding) box for the triangle.
Definition: triangleI.H:213
Point vecC() const
Edge vector opposite point c(): from a() to b()
Definition: triangleI.H:232
Point circumCentre() const
Return circum-centre.
Definition: triangleI.H:287
vector vecA() const
Edge vector opposite point a(): from b() to c()
Definition: triangleI.H:238
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform distribution.
Definition: triangleI.H:406
label n
components
Component labeling enumeration.
Definition: Vector.H:83
volScalarField & p
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:630
Tensor of scalars, i.e. Tensor<scalar>.
const point & b() const noexcept
The second vertex.
Definition: triangle.H:135
vector unitNormal() const
Return unit normal.
Definition: triangleI.H:184
void operator()(const triPoints &)
Definition: triangleI.H:1039
const volScalarField & p0
Definition: EEqn.H:36
point centre() const
Return centre (centroid)
Definition: triangleI.H:132
Namespace for OpenFOAM.
triIntersectionList & tris_
Definition: triangle.H:278
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition: triangleI.H:151
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127