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-2024 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  //- Copy construct from subset of points
114  inline triPoints
115  (
116  const UList<point>& points,
117  const label p0,
118  const label p1,
119  const label p2
120  );
121 
122 
123  // Member Functions
124 
125  //- The first vertex
126  const point& a() const noexcept { return get<0>(); }
127 
128  //- The second vertex
129  const point& b() const noexcept { return get<1>(); }
130 
131  //- The third vertex
132  const point& c() const noexcept { return get<2>(); }
133 
134  //- The first vertex
135  point& a() noexcept { return get<0>(); }
136 
137  //- The second vertex
138  point& b() noexcept { return get<1>(); }
139 
140  //- The third vertex
141  point& c() noexcept { return get<2>(); }
142 
143  //- Flip triangle orientation by swapping second and third vertices
144  inline void flip();
145 
146  //- Return as triangle reference
147  inline triPointRef tri() const;
148 
149 
150  // Properties
152  //- Return centre (centroid)
153  inline point centre() const;
154 
155  //- The area normal - with magnitude equal to area of triangle
156  inline vector areaNormal() const;
157 
158  //- Return unit normal
159  inline vector unitNormal() const;
160 
161  //- The magnitude of the triangle area
162  inline scalar mag() const;
163 
164  //- The enclosing (bounding) box for the triangle
165  inline Pair<point> box() const;
167  //- Edge vector opposite point a(): from b() to c()
168  inline vector vecA() const;
169 
170  //- Edge vector opposite point b(): from c() to a()
171  inline vector vecB() const;
172 
173  //- Edge vector opposite point c(): from a() to b()
174  inline vector vecC() const;
175 };
176 
177 
178 /*---------------------------------------------------------------------------*\
179  Class triangle Declaration
180 \*---------------------------------------------------------------------------*/
181 
182 template<class Point, class PointRef>
183 class triangle
184 {
185 public:
186 
187  // Public Typedefs
188 
189  //- The point type
190  typedef Point point_type;
191 
192  //- Storage type for triangles originating from intersecting triangle
193  //- with another triangle
195 
196  //- Proximity classifications
197  enum proxType
198  {
199  NONE = 0,
200  POINT,
201  EDGE
202  };
203 
204 
205  // Public Classes
206 
207  // Classes for use in sliceWithPlane.
208  // What to do with decomposition of triangle.
209 
210  //- Dummy
211  class dummyOp
212  {
213  public:
214  inline void operator()(const triPoints&);
215  };
216 
217  //- Sum resulting areas
218  class sumAreaOp
219  {
220  public:
221  scalar area_;
222 
223  inline sumAreaOp();
224 
225  inline void operator()(const triPoints&);
226  };
227 
228  //- Store resulting tris
229  class storeOp
230  {
231  public:
233  label& nTris_;
234 
235  inline storeOp(triIntersectionList&, label&);
236 
237  inline void operator()(const triPoints&);
238  };
239 
240 
241 private:
242 
243  // Private Data
244 
245  //- Reference to the first triangle point
246  PointRef a_;
247 
248  //- Reference to the second triangle point
249  PointRef b_;
251  //- Reference to the third triangle point
252  PointRef c_;
253 
254 
255  // Private Member Functions
256 
257  //- Helper: calculate intersection point
258  inline static point planeIntersection
259  (
260  const FixedList<scalar, 3>& d,
261  const triPoints& t,
262  const label negI,
263  const label posI
264  );
265 
266  //- Helper: slice triangle with plane
267  template<class AboveOp, class BelowOp>
268  inline static void triSliceWithPlane
269  (
270  const plane& pln,
271  const triPoints& tri,
272  AboveOp& aboveOp,
273  BelowOp& belowOp
274  );
275 
277 public:
278 
279  // Constructors
280 
281  //- Construct from three points
282  inline triangle(const Point& p0, const Point& p1, const Point& p2);
283 
284  //- Construct from three points
285  inline triangle(const FixedList<Point, 3>& pts);
287  //- Construct from three points out of the list of points
288  // The indices could be from triFace etc.
289  inline triangle
290  (
291  const UList<Point>& points,
292  const FixedList<label, 3>& indices
293  );
294 
295  //- Construct from three points out of the list of points
296  inline triangle
297  (
298  const UList<Point>& points,
299  const label p0,
300  const label p1,
301  const label p2
302  );
303 
304  //- Construct from Istream
305  inline explicit triangle(Istream& is);
306 
307 
308  // Member Functions
309 
310  // Access
311 
312  //- The first vertex
313  const Point& a() const noexcept { return a_; }
314 
315  //- The second vertex
316  const Point& b() const noexcept { return b_; }
317 
318  //- The third vertex
319  const Point& c() const noexcept { return c_; }
320 
321 
322  // Triangle properties (static calculations)
323 
324  //- The centre (centroid) of three points
325  inline static Point centre
326  (
327  const Point& p0,
328  const Point& p1,
329  const Point& p2
330  );
331 
332  //- The area normal for a triangle defined by three points
333  //- (right-hand rule). Magnitude equal to area of the triangle
334  inline static vector areaNormal
335  (
336  const Point& p0,
337  const Point& p1,
338  const Point& p2
339  );
340 
341  //- The unit normal for a triangle defined by three points
342  //- (right-hand rule).
343  inline static vector unitNormal
344  (
345  const Point& p0,
346  const Point& p1,
347  const Point& p2
348  );
349 
350  //- The enclosing (bounding) box for three points
351  inline static Pair<Point> box
352  (
353  const Point& p0,
354  const Point& p1,
355  const Point& p2
356  );
357 
358 
359  // Properties
360 
361  //- Return centre (centroid)
362  inline Point centre() const;
363 
364  //- The area normal - with magnitude equal to area of triangle
365  inline vector areaNormal() const;
366 
367  //- Return unit normal
368  inline vector unitNormal() const;
369 
370  //- Legacy name for areaNormal().
371  // \deprecated(2018-06) Deprecated for new use
372  FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
373  vector normal() const
374  {
375  return areaNormal();
376  }
377 
378  //- The magnitude of the triangle area
379  inline scalar mag() const;
380 
381  //- The enclosing (bounding) box for the triangle
382  inline Pair<Point> box() const;
383 
384  //- Edge vector opposite point a(): from b() to c()
385  inline Point vecA() const;
386 
387  //- Edge vector opposite point b(): from c() to a()
388  inline Point vecB() const;
389 
390  //- Edge vector opposite point c(): from a() to b()
391  inline Point vecC() const;
392 
394  // Properties
395 
396  //- Return circum-centre
397  inline Point circumCentre() const;
399  //- Return circum-radius
400  inline scalar circumRadius() const;
401 
402  //- Return quality: Ratio of triangle and circum-circle
403  // area, scaled so that an equilateral triangle has a
404  // quality of 1
405  inline scalar quality() const;
406 
407  //- Return swept-volume
408  inline scalar sweptVol(const triangle& t) const;
409 
410  //- Return the inertia tensor, with optional reference
411  // point and density specification
412  inline tensor inertia
413  (
414  PointRef refPt = Zero,
415  scalar density = 1.0
416  ) const;
417 
418  //- Return a random point on the triangle from a uniform
419  //- distribution
420  inline Point randomPoint(Random& rndGen) const;
421 
422  //- Calculate the point from the given barycentric coordinates.
423  inline Point barycentricToPoint(const barycentric2D& bary) const;
424 
425  //- Calculate the barycentric coordinates from the given point
426  inline barycentric2D pointToBarycentric(const point& pt) const;
427 
428  //- Calculate the barycentric coordinates from the given point.
429  // Returns the determinant.
430  inline scalar pointToBarycentric
431  (
432  const point& pt,
433  barycentric2D& bary
434  ) const;
435 
436  //- Fast intersection detection with a plane.
437  inline bool intersects
438  (
439  const point& origin,
440  const vector& normal
441  ) const;
442 
443  //- Fast intersection detection with an \b axis plane.
444  inline bool intersects
445  (
447  const point& origin,
449  const vector::components axis
450  ) const;
451 
452  //- Return point intersection with a ray.
453  // For a hit, the distance is signed. Positive number
454  // represents the point in front of triangle.
455  // In case of miss pointHit is set to nearest point
456  // on triangle and its distance to the distance between
457  // the original point and the plane intersection point
458  inline pointHit ray
459  (
460  const point& p,
461  const vector& q,
464  ) const;
465 
466  //- Fast intersection with a ray.
467  // For a hit, the pointHit.distance() is the line parameter t :
468  // intersection=p+t*q. Only defined for VISIBLE, FULL_RAY or
469  // HALF_RAY. tol increases the virtual size of the triangle
470  // by a relative factor.
471  inline pointHit intersection
472  (
473  const point& p,
474  const vector& q,
475  const intersection::algorithm alg,
476  const scalar tol = 0.0
477  ) const;
478 
479  //- Find the nearest point to p on the triangle and classify it:
480  // + near point (nearType=POINT, nearLabel=0, 1, 2)
481  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
482  // Note: edges are counted from starting
483  // vertex so e.g. edge 2 is from f[2] to f[0]
485  (
486  const point& p,
487  label& nearType,
488  label& nearLabel
489  ) const;
490 
491  //- Return nearest point to p on triangle
492  inline pointHit nearestPoint(const point& p) const;
493 
494  //- Classify nearest point to p in triangle plane
495  // w.r.t. triangle edges and points. Returns inside
496  // (true)/outside (false).
497  bool classify
498  (
499  const point& p,
500  label& nearType,
501  label& nearLabel
502  ) const;
503 
504  //- Return nearest point to line on triangle. Returns hit if
505  // point is inside triangle. Sets edgePoint to point on edge
506  // (hit if nearest is inside line)
507  inline pointHit nearestPoint
508  (
509  const linePointRef& edge,
510  pointHit& edgePoint
511  ) const;
512 
513  //- The sign for which side of the face plane the point is on.
514  // Uses the supplied tolerance for rounding around zero.
515  // \return
516  // - 0: on plane
517  // - +1: above plane
518  // - -1: below plane
519  inline int sign(const point& p, const scalar tol = SMALL) const;
520 
521  //- Decompose triangle into triangles above and below plane
522  template<class AboveOp, class BelowOp>
523  inline void sliceWithPlane
524  (
525  const plane& pln,
526  AboveOp& aboveOp,
527  BelowOp& belowOp
528  ) const;
529 
530  //- Decompose triangle into triangles inside and outside
531  // (with respect to user provided normal) other
532  // triangle.
533  template<class InsideOp, class OutsideOp>
534  inline void triangleOverlap
535  (
536  const vector& n,
537  const triangle<Point, PointRef>& tri,
538  InsideOp& insideOp,
539  OutsideOp& outsideOp
540  ) const;
541 
542 
543  // IOstream Operators
544 
545  friend Istream& operator>> <Point, PointRef>
546  (
547  Istream&,
548  triangle&
549  );
550 
551  friend Ostream& operator<< <Point, PointRef>
552  (
553  Ostream&,
554  const triangle&
555  );
556 };
557 
558 
559 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
560 
561 } // End namespace Foam
562 
563 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
564 
565 #include "triangleI.H"
566 
567 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
568 
569 #ifdef NoRepository
570 # include "triangle.C"
571 #endif
572 
573 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
574 
575 #endif
576 
577 // ************************************************************************* //
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition: triangleI.H:380
const point & c() const noexcept
The third vertex.
Definition: triangle.H:151
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:228
A line primitive.
Definition: line.H:52
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: triangleI.H:453
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:309
Unknown proximity.
Definition: triangle.H:250
vector vecC() const
Edge vector opposite point c(): from a() to b()
Definition: triangleI.H:279
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition: triangleI.H:364
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:236
vector vecB() const
Edge vector opposite point b(): from c() to a()
Definition: triangleI.H:273
tensor inertia(PointRef refPt=Zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:399
Random rndGen
Definition: createFields.H:23
const Point & a() const noexcept
The first vertex.
Definition: triangle.H:393
Point point_type
The point type.
Definition: triangle.H:237
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:741
triPointRef tri() const
Return as triangle reference.
Definition: triangleI.H:287
void flip()
Flip triangle orientation by swapping second and third vertices.
Definition: triangleI.H:293
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:403
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:1040
Point vecB() const
Edge vector opposite point b(): from c() to a()
Definition: triangleI.H:255
Close to point.
Definition: triangle.H:251
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:473
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:81
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:904
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:890
proxType
Proximity classifications.
Definition: triangle.H:248
vector unitNormal() const
Return unit normal.
Definition: triangleI.H:207
Close to edge.
Definition: triangle.H:252
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:544
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:398
Point vecA() const
Edge vector opposite point a(): from b() to c()
Definition: triangleI.H:249
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:243
void operator()(const triPoints &)
Definition: triangleI.H:1090
Point centre() const
Return centre (centroid)
Definition: triangleI.H:155
void operator()(const triPoints &)
Definition: triangleI.H:1053
scalar mag() const
The magnitude of the triangle area.
Definition: triangleI.H:302
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:343
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:1078
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:186
const point & a() const noexcept
The first vertex.
Definition: triangle.H:141
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition: triangleI.H:443
bool intersects(const point &origin, const vector &normal) const
Fast intersection detection with a plane.
Definition: triangleI.H:504
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:242
Point vecC() const
Edge vector opposite point c(): from a() to b()
Definition: triangleI.H:261
Point circumCentre() const
Return circum-centre.
Definition: triangleI.H:316
Store resulting tris.
Definition: triangle.H:286
vector vecA() const
Edge vector opposite point a(): from b() to c()
Definition: triangleI.H:267
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform distribution.
Definition: triangleI.H:435
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:659
Tensor of scalars, i.e. Tensor<scalar>.
const point & b() const noexcept
The second vertex.
Definition: triangle.H:146
vector unitNormal() const
Return unit normal.
Definition: triangleI.H:213
void operator()(const triPoints &)
Definition: triangleI.H:1068
const volScalarField & p0
Definition: EEqn.H:36
point centre() const
Return centre (centroid)
Definition: triangleI.H:161
Namespace for OpenFOAM.
triIntersectionList & tris_
Definition: triangle.H:289
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition: triangleI.H:180
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127