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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::triangle
26 
27 Description
28  A triangle primitive used to calculate face normals and swept volumes.
29 
30 SourceFiles
31  triangleI.H
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef triangle_H
36 #define triangle_H
37 
38 #include "intersection.H"
39 #include "vector.H"
40 #include "tensor.H"
41 #include "pointHit.H"
42 #include "Random.H"
43 #include "cachedRandom.H"
44 #include "FixedList.H"
45 #include "UList.H"
46 #include "linePointRef.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 class Istream;
54 class Ostream;
55 class triPoints;
56 class plane;
57 
58 // Forward declaration of friend functions and operators
59 
60 template<class Point, class PointRef> class triangle;
61 
62 template<class Point, class PointRef>
63 inline Istream& operator>>
64 (
65  Istream&,
67 );
68 
69 template<class Point, class PointRef>
70 inline Ostream& operator<<
71 (
72  Ostream&,
74 );
75 
77 
78 
79 /*---------------------------------------------------------------------------*\
80  Class triangle Declaration
81 \*---------------------------------------------------------------------------*/
82 
83 template<class Point, class PointRef>
84 class triangle
85 {
86 public:
87 
88  // Public typedefs
89 
90  //- Storage type for triangles originating from intersecting triangle
91  // with another triangle
93 
94  //- Return types for classify
95  enum proxType
96  {
98  POINT, // Close to point
99  EDGE // Close to edge
100  };
101 
102 
103  // Public classes
104 
105  // Classes for use in sliceWithPlane. What to do with decomposition
106  // of triangle.
107 
108  //- Dummy
109  class dummyOp
110  {
111  public:
112  inline void operator()(const triPoints&);
113  };
114 
115  //- Sum resulting areas
116  class sumAreaOp
117  {
118  public:
119  scalar area_;
120 
121  inline sumAreaOp();
122 
123  inline void operator()(const triPoints&);
124  };
125 
126  //- Store resulting tris
127  class storeOp
128  {
129  public:
131  label& nTris_;
132 
133  inline storeOp(triIntersectionList&, label&);
134 
135  inline void operator()(const triPoints&);
136  };
137 
138 
139 private:
140 
141  // Private data
142 
143  PointRef a_, b_, c_;
144 
145 
146  // Private Member Functions
147 
148  //- Helper: calculate intersection point
149  inline static point planeIntersection
150  (
151  const FixedList<scalar, 3>& d,
152  const triPoints& t,
153  const label negI,
154  const label posI
155  );
156 
157  //- Helper: slice triangle with plane
158  template<class AboveOp, class BelowOp>
159  inline static void triSliceWithPlane
160  (
161  const plane& pl,
162  const triPoints& tri,
163  AboveOp& aboveOp,
164  BelowOp& belowOp
165  );
166 
167 
168 public:
169 
170  // Constructors
171 
172  //- Construct from three points
173  inline triangle(const Point& a, const Point& b, const Point& c);
174 
175  //- Construct from three points
176  inline triangle(const FixedList<Point, 3>&);
177 
178  //- Construct from three points in the list of points
179  // The indices could be from triFace etc.
180  inline triangle
181  (
182  const UList<Point>&,
183  const FixedList<label, 3>& indices
184  );
185 
186  //- Construct from Istream
187  inline triangle(Istream&);
188 
189 
190  // Member Functions
191 
192  // Access
193 
194  //- Return first vertex
195  inline const Point& a() const;
196 
197  //- Return second vertex
198  inline const Point& b() const;
199 
200  //- Return third vertex
201  inline const Point& c() const;
202 
203 
204  // Properties
205 
206  //- Return centre (centroid)
207  inline Point centre() const;
208 
209  //- Return scalar magnitude
210  inline scalar mag() const;
211 
212  //- Return vector normal
213  inline vector normal() const;
214 
215  //- Return circum-centre
216  inline Point circumCentre() const;
217 
218  //- Return circum-radius
219  inline scalar circumRadius() const;
220 
221  //- Return quality: Ratio of triangle and circum-circle
222  // area, scaled so that an equilateral triangle has a
223  // quality of 1
224  inline scalar quality() const;
225 
226  //- Return swept-volume
227  inline scalar sweptVol(const triangle& t) const;
228 
229  //- Return the inertia tensor, with optional reference
230  // point and density specification
231  inline tensor inertia
232  (
233  PointRef refPt = vector::zero,
234  scalar density = 1.0
235  ) const;
236 
237  //- Return a random point on the triangle from a uniform
238  // distribution
239  inline Point randomPoint(Random& rndGen) const;
240 
241  //- Return a random point on the triangle from a uniform
242  // distribution
243  inline Point randomPoint(cachedRandom& rndGen) const;
244 
245  //- Calculate the barycentric coordinates of the given
246  // point, in the same order as a, b, c. Returns the
247  // determinant of the solution.
248  inline scalar barycentric
249  (
250  const point& pt,
251  List<scalar>& bary
252  ) const;
253 
254  //- Return point intersection with a ray.
255  // For a hit, the distance is signed. Positive number
256  // represents the point in front of triangle.
257  // In case of miss pointHit is set to nearest point
258  // on triangle and its distance to the distance between
259  // the original point and the plane intersection point
260  inline pointHit ray
261  (
262  const point& p,
263  const vector& q,
266  ) const;
267 
268  //- Fast intersection with a ray.
269  // For a hit, the pointHit.distance() is the line parameter t :
270  // intersection=p+t*q. Only defined for VISIBLE, FULL_RAY or
271  // HALF_RAY. tol increases the virtual size of the triangle
272  // by a relative factor.
273  inline pointHit intersection
274  (
275  const point& p,
276  const vector& q,
277  const intersection::algorithm alg,
278  const scalar tol = 0.0
279  ) const;
280 
281  //- Find the nearest point to p on the triangle and classify it:
282  // + near point (nearType=POINT, nearLabel=0, 1, 2)
283  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
284  // Note: edges are counted from starting
285  // vertex so e.g. edge 2 is from f[2] to f[0]
287  (
288  const point& p,
289  label& nearType,
290  label& nearLabel
291  ) const;
292 
293  //- Return nearest point to p on triangle
294  inline pointHit nearestPoint(const point& p) const;
295 
296  //- Classify nearest point to p in triangle plane
297  // w.r.t. triangle edges and points. Returns inside
298  // (true)/outside (false).
299  bool classify
300  (
301  const point& p,
302  label& nearType,
303  label& nearLabel
304  ) const;
305 
306  //- Return nearest point to line on triangle. Returns hit if
307  // point is inside triangle. Sets edgePoint to point on edge
308  // (hit if nearest is inside line)
309  inline pointHit nearestPoint
310  (
311  const linePointRef& edge,
312  pointHit& edgePoint
313  ) const;
314 
315  //- Decompose triangle into triangles above and below plane
316  template<class AboveOp, class BelowOp>
317  inline void sliceWithPlane
318  (
319  const plane& pl,
320  AboveOp& aboveOp,
321  BelowOp& belowOp
322  ) const;
323 
324  //- Decompose triangle into triangles inside and outside
325  // (with respect to user provided normal) other
326  // triangle.
327  template<class InsideOp, class OutsideOp>
328  inline void triangleOverlap
329  (
330  const vector& n,
331  const triangle<Point, PointRef>& tri,
332  InsideOp& insideOp,
333  OutsideOp& outsideOp
334  ) const;
335 
336 
337  // IOstream operators
338 
339  friend Istream& operator>> <Point, PointRef>
340  (
341  Istream&,
342  triangle&
343  );
344 
345  friend Ostream& operator<< <Point, PointRef>
346  (
347  Ostream&,
348  const triangle&
349  );
350 };
351 
352 
353 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354 
355 } // End namespace Foam
356 
357 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
358 
359 #include "triangleI.H"
360 
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
362 
363 #ifdef NoRepository
364 # include "triangle.C"
365 #endif
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 #endif
370 
371 // ************************************************************************* //
Foam::intersection::direction
direction
Definition: intersection.H:63
Foam::Random
Simple random number generator.
Definition: Random.H:49
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::triangle::triIntersectionList
FixedList< triPoints, 27 > triIntersectionList
Storage type for triangles originating from intersecting triangle.
Definition: triangle.H:91
Foam::triangle::storeOp::operator()
void operator()(const triPoints &)
Definition: triangleI.H:859
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
p
p
Definition: pEqn.H:62
intersection.H
Foam::triangle::storeOp
Store resulting tris.
Definition: triangle.H:126
Foam::triangle::randomPoint
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform.
Definition: triangleI.H:245
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:51
Foam::triangle::planeIntersection
static point planeIntersection(const FixedList< scalar, 3 > &d, const triPoints &t, const label negI, const label posI)
Helper: calculate intersection point.
Definition: triangleI.H:869
Foam::triangle::sumAreaOp::sumAreaOp
sumAreaOp()
Definition: triangleI.H:829
Foam::triangle::classify
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:684
triangle.C
Foam::triangle::dummyOp
Dummy.
Definition: triangle.H:108
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::intersection::FULL_RAY
@ FULL_RAY
Definition: intersection.H:71
Foam::triangle::a_
PointRef a_
Definition: triangle.H:142
Foam::triangle::sumAreaOp::operator()
void operator()(const triPoints &)
Definition: triangleI.H:837
Foam::triangle::inertia
tensor inertia(PointRef refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:209
Foam::triangle::c_
PointRef c_
Definition: triangle.H:142
pointHit.H
tensor.H
Foam::triangle::centre
Point centre() const
Return centre (centroid)
Definition: triangleI.H:102
Foam::plane
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
Foam::triangle::barycentric
scalar barycentric(const point &pt, List< scalar > &bary) const
Calculate the barycentric coordinates of the given.
Definition: triangleI.H:281
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triPoints
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:49
Foam::triangle::c
const Point & c() const
Return third vertex.
Definition: triangleI.H:95
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:59
Foam::triangle::normal
vector normal() const
Return vector normal.
Definition: triangleI.H:116
Foam::cachedRandom
Random number generator.
Definition: cachedRandom.H:63
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::triangle::a
const Point & a() const
Return first vertex.
Definition: triangleI.H:83
Foam::triangle::NONE
@ NONE
Definition: triangle.H:96
Foam::intersection::algorithm
algorithm
Definition: intersection.H:69
Foam::triangle::triSliceWithPlane
static void triSliceWithPlane(const plane &pl, const triPoints &tri, AboveOp &aboveOp, BelowOp &belowOp)
Helper: slice triangle with plane.
Definition: triangle.C:35
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:322
cachedRandom.H
Foam::triangle::sumAreaOp::area_
scalar area_
Definition: triangle.H:118
Foam::triangle::b
const Point & b() const
Return second vertex.
Definition: triangleI.H:89
Foam::triangle::circumRadius
scalar circumRadius() const
Return circum-radius.
Definition: triangleI.H:150
Foam::triangle::b_
PointRef b_
Definition: triangle.H:142
Foam::triangle::POINT
@ POINT
Definition: triangle.H:97
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::intersection::VECTOR
@ VECTOR
Definition: intersection.H:65
triangleI.H
Random.H
Foam::triangle::sumAreaOp
Sum resulting areas.
Definition: triangle.H:115
Foam::triangle::mag
scalar mag() const
Return scalar magnitude.
Definition: triangleI.H:109
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:521
Foam::triangle::EDGE
@ EDGE
Definition: triangle.H:98
UList.H
Foam::triangle< Foam::Vector, Foam::Vector >::proxType
proxType
Return types for classify.
Definition: triangle.H:94
Foam::Vector< scalar >
Foam::triangle::storeOp::nTris_
label & nTris_
Definition: triangle.H:130
Foam::triangle::quality
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition: triangleI.H:174
Foam::triangle::circumCentre
Point circumCentre() const
Return circum-centre.
Definition: triangleI.H:123
Foam::List< scalar >
Foam::triangle::sliceWithPlane
void sliceWithPlane(const plane &pl, AboveOp &aboveOp, BelowOp &belowOp) const
Decompose triangle into triangles above and below plane.
Definition: triangle.C:113
Foam::triangle::dummyOp::operator()
void operator()(const triPoints &)
Definition: triangleI.H:822
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
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:60
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:438
Foam::line
A line primitive.
Definition: line.H:56
vector.H
Foam::triangle::triangle
triangle(const Point &a, const Point &b, const Point &c)
Construct from three points.
Definition: triangleI.H:35
Foam::triangle::storeOp::tris_
triIntersectionList & tris_
Definition: triangle.H:129
Point
CGAL::Point_3< K > Point
Definition: CGALIndexedPolyhedron.H:51
Foam::triangle::nearestPoint
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:670
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
FixedList.H
rndGen
cachedRandom rndGen(label(0), -1)
Foam::triangle::sweptVol
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition: triangleI.H:190
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:75
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:126
Foam::triangle::storeOp::storeOp
storeOp(triIntersectionList &, label &)
Definition: triangleI.H:847