Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::triangle
29
30Description
31 A triangle primitive used to calculate face normals and swept volumes.
32 Uses referred points.
33
34SourceFiles
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
59namespace Foam
60{
62// Forward Declarations
63class plane;
65template<class Point, class PointRef>
67
68template<class Point, class PointRef>
70
71
72/*---------------------------------------------------------------------------*\
73 Class triPoints Declaration
74\*---------------------------------------------------------------------------*/
75
76//- Triangle point storage. Default constructable (triangle is not)
77class triPoints
78:
79 public FixedList<point, 3>
80{
81public:
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;
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 magnitude squared of the triangle area
165 inline scalar magSqr() const;
167 //- The enclosing (bounding) box for the triangle
168 inline Pair<point> box() const;
169
170 //- Edge vector opposite point a(): from b() to c()
171 inline vector vecA() const;
172
173 //- Edge vector opposite point b(): from c() to a()
174 inline vector vecB() const;
175
176 //- Edge vector opposite point c(): from a() to b()
177 inline vector vecC() const;
178};
179
180
181/*---------------------------------------------------------------------------*\
182 Class triangle Declaration
183\*---------------------------------------------------------------------------*/
184
185template<class Point, class PointRef>
186class triangle
187{
188public:
189
190 // Public Typedefs
191
192 //- The point type
193 typedef Point point_type;
194
195 //- Storage type for triangles originating from intersecting triangle
196 //- with another triangle
198
199 //- Proximity classifications
200 enum proxType
201 {
202 NONE = 0,
203 POINT,
204 EDGE
205 };
206
207
208 // Public Classes
209
210 // Classes for use in sliceWithPlane.
211 // What to do with decomposition of triangle.
212
213 //- Dummy
214 struct dummyOp
215 {
216 void operator()(const triPoints&) const noexcept
217 {}
218 };
219
220 //- Sum resulting areas
221 struct sumAreaOp
222 {
223 scalar area_;
224
225 constexpr sumAreaOp() noexcept : area_(0) {}
226
227 scalar total() const noexcept { return area_; }
228
229 inline void operator()(const triPoints&);
230 };
231
232 //- Store resulting tris
233 struct storeOp
234 {
236 label& nTris_;
237
238 inline storeOp(triIntersectionList&, label&);
239
240 inline void operator()(const triPoints&);
241 };
243
244private:
245
246 // Private Data
247
248 //- Reference to the first triangle point
249 PointRef a_;
250
251 //- Reference to the second triangle point
252 PointRef b_;
254 //- Reference to the third triangle point
255 PointRef c_;
258 // Private Member Functions
259
260 //- Helper: calculate intersection point
261 inline static point planeIntersection
262 (
263 const FixedList<scalar, 3>& d,
264 const triPoints& t,
265 const label negI,
266 const label posI
267 );
268
269 //- Helper: slice triangle with plane
270 template<class AboveOp, class BelowOp>
271 inline static void triSliceWithPlane
272 (
273 const plane& pln,
274 const triPoints& tri,
275 AboveOp& aboveOp,
276 BelowOp& belowOp
277 );
279
280public:
281
282 // Constructors
283
284 //- Construct from three points
285 inline triangle(const Point& p0, const Point& p1, const Point& p2);
286
287 //- Construct from three points
288 inline triangle(const FixedList<Point, 3>& pts);
289
290 //- Construct from three points out of the list of points
291 // The indices could be from triFace etc.
292 inline triangle
293 (
295 const FixedList<label, 3>& indices
296 );
297
298 //- Construct from three points out of the list of points
299 inline triangle
300 (
301 const UList<Point>& points,
302 const label p0,
303 const label p1,
304 const label p2
305 );
306
307 //- Construct from Istream
308 inline explicit triangle(Istream& is);
309
310
311 // Member Functions
312
313 // Access
314
315 //- The first vertex
316 const Point& a() const noexcept { return a_; }
317
318 //- The second vertex
319 const Point& b() const noexcept { return b_; }
320
321 //- The third vertex
322 const Point& c() const noexcept { return c_; }
323
324
325 // Triangle properties (static calculations)
326
327 //- The centre (centroid) of three points
328 inline static Point centre
329 (
330 const Point& p0,
331 const Point& p1,
332 const Point& p2
333 );
334
335 //- The area normal for a triangle defined by three points
336 //- (right-hand rule). Magnitude equal to area of the triangle
337 inline static vector areaNormal
338 (
339 const Point& p0,
340 const Point& p1,
341 const Point& p2
342 );
343
344 //- The unit normal for a triangle defined by three points
345 //- (right-hand rule).
346 inline static vector unitNormal
347 (
348 const Point& p0,
349 const Point& p1,
350 const Point& p2
351 );
352
353 //- The enclosing (bounding) box for three points
354 inline static Pair<Point> box
355 (
356 const Point& p0,
357 const Point& p1,
358 const Point& p2
359 );
360
362 // Properties
363
364 //- Return centre (centroid)
365 inline Point centre() const;
366
367 //- The area normal - with magnitude equal to area of triangle
368 inline vector areaNormal() const;
369
370 //- Return unit normal
371 inline vector unitNormal() const;
372
373 //- Legacy name for areaNormal().
374 // \deprecated(2018-06) Deprecated for new use
375 FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()")
376 vector normal() const
378 return areaNormal();
379 }
380
381 //- The magnitude of the triangle area
382 inline scalar mag() const;
383
384 //- The magnitude squared of the triangle area
385 inline scalar magSqr() const;
386
387 //- The enclosing (bounding) box for the triangle
388 inline Pair<Point> box() const;
389
390 //- Edge vector opposite point a(): from b() to c()
391 inline Point vecA() const;
392
393 //- Edge vector opposite point b(): from c() to a()
394 inline Point vecB() const;
395
396 //- Edge vector opposite point c(): from a() to b()
397 inline Point vecC() const;
399
400 // Properties
401
402 //- Return circum-centre
403 inline Point circumCentre() const;
404
405 //- Return circum-radius
406 inline scalar circumRadius() const;
407
408 //- Return quality: Ratio of triangle and circum-circle
409 // area, scaled so that an equilateral triangle has a
410 // quality of 1
411 inline scalar quality() const;
412
413 //- Return swept-volume
414 inline scalar sweptVol(const triangle& t) const;
415
416 //- Return the inertia tensor, with optional reference
417 // point and density specification
418 inline tensor inertia
419 (
420 PointRef refPt = Zero,
421 scalar density = 1.0
422 ) const;
423
424 //- Return a random point on the triangle from a uniform
425 //- distribution
426 inline Point randomPoint(Random& rndGen) const;
428 //- Calculate the point from the given barycentric coordinates.
429 inline Point barycentricToPoint(const barycentric2D& bary) const;
430
431 //- Calculate the barycentric coordinates from the given point
432 inline barycentric2D pointToBarycentric(const point& pt) const;
433
434 //- Calculate the barycentric coordinates from the given point.
435 // Returns the determinant.
436 inline scalar pointToBarycentric
437 (
438 const point& pt,
439 barycentric2D& bary
440 ) const;
441
442 //- Fast intersection detection with a plane.
443 inline bool intersects
444 (
445 const point& origin,
446 const vector& normal
447 ) const;
449 //- Fast intersection detection with an \b axis plane.
450 inline bool intersects
451 (
453 const point& origin,
455 const vector::components axis
456 ) const;
457
458 //- Return point intersection with a ray.
459 // For a hit, the distance is signed. Positive number
460 // represents the point in front of triangle.
461 // In case of miss pointHit is set to nearest point
462 // on triangle and its distance to the distance between
463 // the original point and the plane intersection point
464 inline pointHit ray
465 (
466 const point& p,
467 const vector& q,
470 ) const;
472 //- Fast intersection with a ray.
473 // For a hit, the pointHit.distance() is the line parameter t :
474 // intersection=p+t*q. Only defined for VISIBLE, FULL_RAY or
475 // HALF_RAY. tol increases the virtual size of the triangle
476 // by a relative factor.
479 const point& p,
480 const vector& q,
481 const intersection::algorithm alg,
482 const scalar tol = 0.0
483 ) const;
484
485 //- Find the nearest point to p on the triangle and classify it:
486 // + near point (nearType=POINT, nearLabel=0, 1, 2)
487 // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
488 // Note: edges are counted from starting
489 // vertex so e.g. edge 2 is from f[2] to f[0]
491 (
492 const point& p,
493 label& nearType,
494 label& nearLabel
495 ) const;
496
497 //- Return nearest point to p on triangle
498 inline pointHit nearestPoint(const point& p) const;
499
500 //- Classify nearest point to p in triangle plane
501 // w.r.t. triangle edges and points. Returns inside
502 // (true)/outside (false).
503 bool classify
504 (
505 const point& p,
506 label& nearType,
507 label& nearLabel
508 ) const;
509
510 //- Return nearest point to line on triangle. Returns hit if
511 // point is inside triangle. Sets edgePoint to point on edge
512 // (hit if nearest is inside line)
514 (
515 const linePointRef& edge,
516 pointHit& edgePoint
517 ) const;
518
519 //- The sign for which side of the face plane the point is on.
520 // Uses the supplied tolerance for rounding around zero.
521 // \return
522 // - 0: on plane
523 // - +1: above plane
524 // - -1: below plane
525 inline int sign(const point& p, const scalar tol = SMALL) const;
526
527 //- Decompose triangle into triangles above and below plane
528 template<class AboveOp, class BelowOp>
529 inline void sliceWithPlane
530 (
531 const plane& pln,
532 AboveOp& aboveOp,
533 BelowOp& belowOp
534 ) const;
535
536 //- Decompose triangle into triangles inside and outside
537 // (with respect to user provided normal) other
538 // triangle.
539 template<class InsideOp, class OutsideOp>
540 inline void triangleOverlap
541 (
542 const vector& n,
543 const triangle<Point, PointRef>& tri,
544 InsideOp& insideOp,
545 OutsideOp& outsideOp
546 ) const;
547
548
549 // IOstream Operators
550
551 friend Istream& operator>> <Point, PointRef>
552 (
553 Istream&,
554 triangle&
555 );
556
558 (
559 Ostream&,
560 const triangle&
561 );
562};
563
564
565// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
566
567} // End namespace Foam
568
569// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
570
571#include "triangleI.H"
573// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
574
575#ifdef NoRepository
576# include "triangle.C"
577#endif
578
579// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
580
581#endif
582
583// ************************************************************************* //
CGAL::Point_3< K > Point
label n
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
point & get() noexcept
Definition FixedListI.H:152
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
Random number generator.
Definition Random.H:56
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
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:62
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
Triangle point storage. Default constructable (triangle is not).
Definition triangle.H:77
void flip()
Flip triangle orientation by swapping second and third vertices.
Definition triangleI.H:293
point & b() noexcept
The second vertex.
Definition triangle.H:161
triPointRef tri() const
Return as triangle reference.
Definition triangleI.H:287
void back()=delete
The back() accessor (from FixedList) has no purpose.
void front()=delete
The front() accessor (from FixedList) has no purpose.
const point & c() const noexcept
The third vertex.
Definition triangle.H:151
point & a() noexcept
The first vertex.
Definition triangle.H:156
vector vecA() const
Edge vector opposite point a(): from b() to c().
Definition triangleI.H:267
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition triangleI.H:186
const point & b() const noexcept
The second vertex.
Definition triangle.H:146
Pair< point > box() const
The enclosing (bounding) box for the triangle.
Definition triangleI.H:242
vector unitNormal() const
Return unit normal.
Definition triangleI.H:213
scalar mag() const
The magnitude of the triangle area.
Definition triangleI.H:302
vector vecC() const
Edge vector opposite point c(): from a() to b().
Definition triangleI.H:279
triPoints()=default
Default construct.
point & c() noexcept
The third vertex.
Definition triangle.H:166
const point & a() const noexcept
The first vertex.
Definition triangle.H:141
scalar magSqr() const
The magnitude squared of the triangle area.
Definition triangleI.H:315
point centre() const
Return centre (centroid).
Definition triangleI.H:161
vector vecB() const
Edge vector opposite point b(): from c() to a().
Definition triangleI.H:273
A triangle primitive used to calculate face normals and swept volumes. Uses referred points.
Definition triangle.H:234
Point circumCentre() const
Return circum-centre.
Definition triangleI.H:329
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition triangleI.H:466
void sliceWithPlane(const plane &pln, AboveOp &aboveOp, BelowOp &belowOp) const
Decompose triangle into triangles above and below plane.
Definition triangle.C:107
static vector areaNormal(const Point &p0, const Point &p1, const Point &p2)
The area normal for a triangle defined by three points (right-hand rule). Magnitude equal to area of ...
Definition triangleI.H:169
triangle(const UList< Point > &points, const label p0, const label p1, const label p2)
Construct from three points out of the list of points.
Definition triangleI.H:120
FOAM_DEPRECATED_FOR(2018-12, "areaNormal() or unitNormal()") vector normal() const
Legacy name for areaNormal().
Definition triangle.H:478
scalar circumRadius() const
Return circum-radius.
Definition triangleI.H:356
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition triangleI.H:456
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition triangleI.H:393
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition triangleI.H:180
Point centre() const
Return centre (centroid).
Definition triangleI.H:155
triangle(const FixedList< Point, 3 > &pts)
Construct from three points.
Definition triangleI.H:95
Point vecB() const
Edge vector opposite point b(): from c() to a().
Definition triangleI.H:255
FixedList< triPoints, 27 > triIntersectionList
Storage type for triangles originating from intersecting triangle with another triangle.
Definition triangle.H:248
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
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:754
bool intersects(const point &origin, const vector::components axis) const
Fast intersection detection with an axis plane.
Definition triangleI.H:537
proxType
Proximity classifications.
Definition triangle.H:254
triangle(const Point &p0, const Point &p1, const Point &p2)
Construct from three points.
Definition triangleI.H:81
static Pair< Point > box(const Point &p0, const Point &p1, const Point &p2)
The enclosing (bounding) box for three points.
Definition triangleI.H:221
static Point centre(const Point &p0, const Point &p1, const Point &p2)
The centre (centroid) of three points.
Definition triangleI.H:144
const Point & b() const noexcept
The second vertex.
Definition triangle.H:403
scalar pointToBarycentric(const point &pt, barycentric2D &bary) const
Calculate the barycentric coordinates from the given point.
Definition triangleI.H:478
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:557
vector unitNormal() const
Return unit normal.
Definition triangleI.H:207
triangle(Istream &is)
Construct from Istream.
Definition triangleI.H:134
pointHit nearestPoint(const linePointRef &edge, pointHit &edgePoint) const
Return nearest point to line on triangle. Returns hit if.
Definition triangleI.H:929
tensor inertia(PointRef refPt=Zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition triangleI.H:412
scalar mag() const
The magnitude of the triangle area.
Definition triangleI.H:309
bool intersects(const point &origin, const vector &normal) const
Fast intersection detection with a plane.
Definition triangleI.H:517
friend Ostream & operator(Ostream &, const triangle &)
Point point_type
The point type.
Definition triangle.H:242
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition triangleI.H:903
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform distribution.
Definition triangleI.H:448
Pair< Point > box() const
The enclosing (bounding) box for the triangle.
Definition triangleI.H:236
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:672
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition triangleI.H:917
scalar magSqr() const
The magnitude squared of the triangle area.
Definition triangleI.H:322
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:1053
static vector unitNormal(const Point &p0, const Point &p1, const Point &p2)
The unit normal for a triangle defined by three points (right-hand rule).
Definition triangleI.H:194
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition triangleI.H:377
const Point & c() const noexcept
The third vertex.
Definition triangle.H:408
Point vecA() const
Edge vector opposite point a(): from b() to c().
Definition triangleI.H:249
const Point & a() const noexcept
The first vertex.
Definition triangle.H:398
triangle(const UList< Point > &points, const FixedList< label, 3 > &indices)
Construct from three points out of the list of points.
Definition triangleI.H:107
Point vecC() const
Edge vector opposite point c(): from a() to b().
Definition triangleI.H:261
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
const pointField & points
Namespace for OpenFOAM.
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
Barycentric2D< scalar > barycentric2D
A scalar version of the templated Barycentric2D.
Tensor< scalar > tensor
Definition symmTensor.H:57
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Istream & operator>>(Istream &, directionInfo &)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const direction noexcept
Definition scalarImpl.H:265
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
Vector< scalar > vector
Definition vector.H:57
const pointField & pts
#define FOAM_DEPRECATED_FOR(since, replacement)
Definition stdFoam.H:43
void operator()(const triPoints &) const noexcept
Definition triangle.H:271
Store resulting tris.
Definition triangle.H:293
storeOp(triIntersectionList &, label &)
Definition triangleI.H:1076
triIntersectionList & tris_
Definition triangle.H:294
void operator()(const triPoints &)
Definition triangleI.H:1088
Sum resulting areas.
Definition triangle.H:279
void operator()(const triPoints &)
Definition triangleI.H:1066
constexpr sumAreaOp() noexcept
Definition triangle.H:282
scalar total() const noexcept
Definition triangle.H:284
Random rndGen