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-2017 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 #include "barycentric2D.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 class Istream;
55 class Ostream;
56 
57 // Forward declaration of friend functions and operators
58 
59 template<class Point, class PointRef> class triangle;
60 
61 template<class Point, class PointRef>
62 inline Istream& operator>>
63 (
64  Istream&,
66 );
67 
68 template<class Point, class PointRef>
69 inline Ostream& operator<<
70 (
71  Ostream&,
73 );
74 
75 
76 /*---------------------------------------------------------------------------*\
77  Class triangle Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 template<class Point, class PointRef>
81 class triangle
82 {
83  // Private data
84 
85  PointRef a_, b_, c_;
86 
87 
88 public:
89 
90  //- Return types for classify
91  enum proxType
92  {
94  POINT, // Close to point
95  EDGE // Close to edge
96  };
97 
98 
99  // Constructors
100 
101  //- Construct from three points
102  inline triangle(const Point& a, const Point& b, const Point& c);
103 
104  //- Construct from three points in the list of points
105  // The indices could be from triFace etc.
106  inline triangle
107  (
108  const UList<Point>&,
109  const FixedList<label, 3>& indices
110  );
111 
112  //- Construct from Istream
113  inline triangle(Istream&);
114 
115 
116  // Member Functions
117 
118  // Access
119 
120  //- Return first vertex
121  inline const Point& a() const;
122 
123  //- Return second vertex
124  inline const Point& b() const;
125 
126  //- Return third vertex
127  inline const Point& c() const;
128 
129 
130  // Properties
131 
132  //- Return centre (centroid)
133  inline Point centre() const;
134 
135  //- Return scalar magnitude
136  inline scalar mag() const;
137 
138  //- Return vector normal
139  inline vector normal() const;
140 
141  //- Return circum-centre
142  inline Point circumCentre() const;
143 
144  //- Return circum-radius
145  inline scalar circumRadius() const;
146 
147  //- Return quality: Ratio of triangle and circum-circle
148  // area, scaled so that an equilateral triangle has a
149  // quality of 1
150  inline scalar quality() const;
151 
152  //- Return swept-volume
153  inline scalar sweptVol(const triangle& t) const;
154 
155  //- Return the inertia tensor, with optional reference
156  // point and density specification
157  inline tensor inertia
158  (
159  PointRef refPt = vector::zero,
160  scalar density = 1.0
161  ) const;
162 
163  //- Return a random point on the triangle from a uniform
164  // distribution
165  inline Point randomPoint(Random& rndGen) const;
166 
167  //- Return a random point on the triangle from a uniform
168  // distribution
169  inline Point randomPoint(cachedRandom& rndGen) const;
170 
171  //- Calculate the point from the given barycentric coordinates.
172  inline Point barycentricToPoint(const barycentric2D& bary) const;
173 
174  //- Calculate the barycentric coordinates from the given point
175  inline barycentric2D pointToBarycentric(const point& pt) const;
176 
177  //- Calculate the barycentric coordinates from the given point.
178  // Returns the determinant.
179  inline scalar pointToBarycentric
180  (
181  const point& pt,
182  barycentric2D& bary
183  ) const;
184 
185  //- Return point intersection with a ray.
186  // For a hit, the distance is signed. Positive number
187  // represents the point in front of triangle.
188  // In case of miss pointHit is set to nearest point
189  // on triangle and its distance to the distance between
190  // the original point and the plane intersection point
191  inline pointHit ray
192  (
193  const point& p,
194  const vector& q,
197  ) const;
198 
199  //- Fast intersection with a ray.
200  // For a hit, the pointHit.distance() is the line parameter t :
201  // intersection=p+t*q. Only defined for VISIBLE, FULL_RAY or
202  // HALF_RAY. tol increases the virtual size of the triangle
203  // by a relative factor.
204  inline pointHit intersection
205  (
206  const point& p,
207  const vector& q,
208  const intersection::algorithm alg,
209  const scalar tol = 0.0
210  ) const;
211 
212  //- Find the nearest point to p on the triangle and classify it:
213  // + near point (nearType=POINT, nearLabel=0, 1, 2)
214  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
215  // Note: edges are counted from starting
216  // vertex so e.g. edge 2 is from f[2] to f[0]
218  (
219  const point& p,
220  label& nearType,
221  label& nearLabel
222  ) const;
223 
224  //- Return nearest point to p on triangle
225  inline pointHit nearestPoint(const point& p) const;
226 
227  //- Classify nearest point to p in triangle plane
228  // w.r.t. triangle edges and points. Returns inside
229  // (true)/outside (false).
230  bool classify
231  (
232  const point& p,
233  label& nearType,
234  label& nearLabel
235  ) const;
236 
237  //- Return nearest point to line on triangle. Returns hit if
238  // point is inside triangle. Sets edgePoint to point on edge
239  // (hit if nearest is inside line)
240  inline pointHit nearestPoint
241  (
242  const linePointRef& edge,
243  pointHit& edgePoint
244  ) const;
245 
246 
247  // IOstream operators
248 
249  friend Istream& operator>> <Point, PointRef>
250  (
251  Istream&,
252  triangle&
253  );
254 
255  friend Ostream& operator<< <Point, PointRef>
256  (
257  Ostream&,
258  const triangle&
259  );
260 };
261 
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 } // End namespace Foam
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 #include "triangleI.H"
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 #endif
274 
275 // ************************************************************************* //
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition: triangleI.H:177
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:58
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
A line primitive.
Definition: line.H:56
const Point & c() const
Return third vertex.
Definition: triangleI.H:82
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: triangleI.H:260
scalar mag() const
Return scalar magnitude.
Definition: triangleI.H:96
scalar quality() const
Return quality: Ratio of triangle and circum-circle.
Definition: triangleI.H:161
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
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:510
Random number generator.
Definition: cachedRandom.H:63
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant...
Definition: Barycentric2D.H:50
bool classify(const point &p, label &nearType, label &nearLabel) const
Classify nearest point to p in triangle plane.
Definition: triangleI.H:673
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:659
proxType
Return types for classify.
Definition: triangle.H:90
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:311
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
cachedRandom rndGen(label(0), -1)
Point centre() const
Return centre (centroid)
Definition: triangleI.H:89
scalar circumRadius() const
Return circum-radius.
Definition: triangleI.H:137
tensor inertia(PointRef refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:196
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:61
Simple random number generator.
Definition: Random.H:49
vector normal() const
Return vector normal.
Definition: triangleI.H:103
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
triangle(const Point &a, const Point &b, const Point &c)
Construct from three points.
Definition: triangleI.H:34
Point barycentricToPoint(const barycentric2D &bary) const
Calculate the point from the given barycentric coordinates.
Definition: triangleI.H:250
const Point & a() const
Return first vertex.
Definition: triangleI.H:70
Point circumCentre() const
Return circum-centre.
Definition: triangleI.H:110
const Point & b() const
Return second vertex.
Definition: triangleI.H:76
Point randomPoint(Random &rndGen) const
Return a random point on the triangle from a uniform.
Definition: triangleI.H:232
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:427
Namespace for OpenFOAM.