tetrahedron.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-2016 OpenFOAM Foundation
9  Copyright (C) 2018 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::tetrahedron
29 
30 Description
31  A tetrahedron primitive.
32 
33  Ordering of edges needs to be the same for a tetrahedron
34  class, a tetrahedron cell shape model and a tetCell.
35 
36 SourceFiles
37  tetrahedronI.H
38  tetrahedron.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef tetrahedron_H
43 #define tetrahedron_H
44 
45 #include "point.H"
46 #include "primitiveFieldsFwd.H"
47 #include "pointHit.H"
48 #include "Random.H"
49 #include "FixedList.H"
50 #include "UList.H"
51 #include "triPointRef.H"
52 #include "boundBox.H"
53 #include "barycentric.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 
60 // Forward declarations
61 class Istream;
62 class Ostream;
63 class tetPoints;
64 class plane;
65 
66 template<class Point, class PointRef> class tetrahedron;
67 
68 template<class Point, class PointRef>
69 inline Istream& operator>>
70 (
71  Istream&,
73 );
74 
75 template<class Point, class PointRef>
76 inline Ostream& operator<<
77 (
78  Ostream&,
80 );
81 
83 
84 /*---------------------------------------------------------------------------*\
85  class tetrahedron Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 template<class Point, class PointRef>
89 class tetrahedron
90 {
91 public:
92 
93  // Public typedefs
94 
95  //- Storage type for tets originating from intersecting tets.
96  // (can possibly be smaller than 200)
98 
99 
100  // Classes for use in sliceWithPlane. What to do with decomposition
101  // of tet.
102 
103  //- Dummy
104  class dummyOp
105  {
106  public:
107  inline void operator()(const tetPoints&);
108  };
109 
110  //- Sum resulting volumes
111  class sumVolOp
112  {
113  public:
114  scalar vol_;
115 
116  inline sumVolOp();
117 
118  inline void operator()(const tetPoints&);
119  };
120 
121  //- Store resulting tets
122  class storeOp
123  {
124  tetIntersectionList& tets_;
125  label& nTets_;
126 
127  public:
128  inline storeOp(tetIntersectionList&, label&);
129 
130  inline void operator()(const tetPoints&);
131  };
132 
133 private:
134 
135  // Private data
136 
137  PointRef a_, b_, c_, d_;
138 
139  inline static point planeIntersection
140  (
141  const FixedList<scalar, 4>&,
142  const tetPoints&,
143  const label,
144  const label
145  );
146 
147  template<class TetOp>
148  inline static void decomposePrism
149  (
151  TetOp& op
152  );
153 
154  template<class AboveTetOp, class BelowTetOp>
155  inline static void tetSliceWithPlane
156  (
157  const plane& pl,
158  const tetPoints& tet,
159  AboveTetOp& aboveOp,
160  BelowTetOp& belowOp
161  );
162 
163 
164 public:
165 
166  // Member constants
167 
168  enum
169  {
170  nVertices = 4, // Number of vertices in tetrahedron
171  nEdges = 6 // Number of edges in tetrahedron
172  };
173 
174 
175  // Constructors
176 
177  //- Construct from points
178  inline tetrahedron
179  (
180  const Point& a,
181  const Point& b,
182  const Point& c,
183  const Point& d
184  );
185 
186  //- Construct from four points in the list of points
187  inline tetrahedron
188  (
189  const UList<Point>&,
190  const FixedList<label, 4>& indices
191  );
192 
193  //- Construct from Istream
194  inline tetrahedron(Istream&);
195 
196 
197  // Member Functions
198 
199  // Access
200 
201  //- Return vertices
202  inline const Point& a() const;
203 
204  inline const Point& b() const;
205 
206  inline const Point& c() const;
207 
208  inline const Point& d() const;
209 
210  //- Return i-th face
211  inline triPointRef tri(const label facei) const;
212 
213  // Properties
214 
215  //- Face area normal for side a
216  inline vector Sa() const;
217 
218  //- Face area normal for side b
219  inline vector Sb() const;
220 
221  //- Face area normal for side c
222  inline vector Sc() const;
223 
224  //- Face area normal for side d
225  inline vector Sd() const;
226 
227  //- Return centre (centroid)
228  inline Point centre() const;
229 
230  //- Return volume
231  inline scalar mag() const;
232 
233  //- Return circum-centre
234  inline Point circumCentre() const;
235 
236  //- Return circum-radius
237  inline scalar circumRadius() const;
238 
239  //- Return quality: Ratio of tetrahedron and circum-sphere
240  // volume, scaled so that a regular tetrahedron has a
241  // quality of 1
242  inline scalar quality() const;
243 
244  //- Return a random point in the tetrahedron from a
245  // uniform distribution
246  inline Point randomPoint(Random& rndGen) const;
247 
248  //- Calculate the point from the given barycentric coordinates.
249  inline Point barycentricToPoint(const barycentric& bary) const;
250 
251  //- Calculate the barycentric coordinates from the given point
252  inline barycentric pointToBarycentric(const point& pt) const;
253 
254  //- Calculate the barycentric coordinates from the given point.
255  // Returns the determinant.
256  inline scalar pointToBarycentric
257  (
258  const point& pt,
259  barycentric& bary
260  ) const;
261 
262  //- Return nearest point to p on tetrahedron. Is p itself
263  // if inside.
264  inline pointHit nearestPoint(const point& p) const;
265 
266  //- Return true if point is inside tetrahedron
267  inline bool inside(const point& pt) const;
268 
269  //- Decompose tet into tets above and below plane
270  template<class AboveTetOp, class BelowTetOp>
271  inline void sliceWithPlane
272  (
273  const plane& pl,
274  AboveTetOp& aboveOp,
275  BelowTetOp& belowOp
276  ) const;
277 
278  //- Decompose tet into tets inside and outside other tet
279  inline void tetOverlap
280  (
281  const tetrahedron<Point, PointRef>& tetB,
282  tetIntersectionList& insideTets,
283  label& nInside,
284  tetIntersectionList& outsideTets,
285  label& nOutside
286  ) const;
287 
288 
289  //- Return (min)containment sphere, i.e. the smallest sphere with
290  // all points inside. Returns pointHit with:
291  // - hit : if sphere is equal to circumsphere
292  // (biggest sphere)
293  // - point : centre of sphere
294  // - distance : radius of sphere
295  // - eligiblemiss: false
296  // Tol (small compared to 1, e.g. 1e-9) is used to determine
297  // whether point is inside: mag(pt - ctr) < (1+tol)*radius.
298  pointHit containmentSphere(const scalar tol) const;
299 
300  //- Fill buffer with shape function products
301  void gradNiSquared(scalarField& buffer) const;
302 
303  void gradNiDotGradNj(scalarField& buffer) const;
304 
305  void gradNiGradNi(tensorField& buffer) const;
306 
307  void gradNiGradNj(tensorField& buffer) const;
308 
309  //- Calculate the bounding box
310  boundBox bounds() const;
311 
312 
313  // IOstream operators
314 
315  friend Istream& operator>> <Point, PointRef>
316  (
317  Istream&,
318  tetrahedron&
319  );
320 
321  friend Ostream& operator<< <Point, PointRef>
322  (
323  Ostream&,
324  const tetrahedron&
325  );
326 };
327 
328 
329 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330 
331 } // End namespace Foam
332 
333 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334 
335 #include "tetrahedronI.H"
336 
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338 
339 #ifdef NoRepository
340  #include "tetrahedron.C"
341 #endif
342 
343 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344 
345 #endif
346 
347 // ************************************************************************* //
Foam::tetrahedron::inside
bool inside(const point &pt) const
Return true if point is inside tetrahedron.
Definition: tetrahedronI.H:415
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::tetrahedron::nearestPoint
pointHit nearestPoint(const point &p) const
Return nearest point to p on tetrahedron. Is p itself.
Definition: tetrahedronI.H:321
primitiveFieldsFwd.H
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
Foam::tetrahedron::storeOp::operator()
void operator()(const tetPoints &)
Definition: tetrahedronI.H:529
Foam::tetrahedron::d
const Point & d() const
Definition: tetrahedronI.H:97
point.H
Foam::tetrahedron::dummyOp
Dummy.
Definition: tetrahedron.H:103
Foam::tetrahedron::sumVolOp
Sum resulting volumes.
Definition: tetrahedron.H:110
Foam::tetrahedron::sumVolOp::sumVolOp
sumVolOp()
Definition: tetrahedronI.H:499
Foam::tetrahedron::circumCentre
Point circumCentre() const
Return circum-centre.
Definition: tetrahedronI.H:179
Foam::tetrahedron::tetrahedron
tetrahedron(const Point &a, const Point &b, const Point &c, const Point &d)
Construct from points.
Definition: tetrahedronI.H:38
triPointRef.H
pointHit.H
Foam::tetrahedron::gradNiGradNi
void gradNiGradNi(tensorField &buffer) const
Definition: tetrahedron.C:299
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetPointRef.H:47
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:61
Foam::tetrahedron::c
const Point & c() const
Definition: tetrahedronI.H:90
Foam::tetrahedron::gradNiDotGradNj
void gradNiDotGradNj(scalarField &buffer) const
Definition: tetrahedron.C:269
Foam::Field< scalar >
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::tetrahedron::containmentSphere
pointHit containmentSphere(const scalar tol) const
Return (min)containment sphere, i.e. the smallest sphere with.
Definition: tetrahedron.C:36
Foam::tetrahedron::Sb
vector Sb() const
Face area normal for side b.
Definition: tetrahedronI.H:144
Foam::Barycentric< scalar >
Foam::tetrahedron::dummyOp::operator()
void operator()(const tetPoints &)
Definition: tetrahedronI.H:492
Foam::tetrahedron::nVertices
Definition: tetrahedron.H:169
Foam::tetrahedron::circumRadius
scalar circumRadius() const
Return circum-radius.
Definition: tetrahedronI.H:206
Foam::tetrahedron::quality
scalar quality() const
Return quality: Ratio of tetrahedron and circum-sphere.
Definition: tetrahedronI.H:232
Foam::tetrahedron::gradNiSquared
void gradNiSquared(scalarField &buffer) const
Fill buffer with shape function products.
Definition: tetrahedron.C:248
Foam::tetrahedron::sumVolOp::vol_
scalar vol_
Definition: tetrahedron.H:113
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::tetrahedron::storeOp
Store resulting tets.
Definition: tetrahedron.H:121
Random.H
Foam::tetrahedron::pointToBarycentric
barycentric pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: tetrahedronI.H:266
Foam::tetPoints
Tet storage. Null constructable (unfortunately tetrahedron<point, point> is not)
Definition: tetPoints.H:53
Foam::tetrahedron::a
const Point & a() const
Return vertices.
Definition: tetrahedronI.H:76
boundBox.H
Foam::tetrahedron::tetOverlap
void tetOverlap(const tetrahedron< Point, PointRef > &tetB, tetIntersectionList &insideTets, label &nInside, tetIntersectionList &outsideTets, label &nOutside) const
Decompose tet into tets inside and outside other tet.
Foam::tetrahedron::gradNiGradNj
void gradNiGradNj(tensorField &buffer) const
Definition: tetrahedron.C:317
Foam::tetrahedron::tri
triPointRef tri(const label facei) const
Return i-th face.
Definition: tetrahedronI.H:105
Foam::tetrahedron::Sa
vector Sa() const
Face area normal for side a.
Definition: tetrahedronI.H:137
UList.H
Foam::tetrahedron::randomPoint
Point randomPoint(Random &rndGen) const
Return a random point in the tetrahedron from a.
Definition: tetrahedronI.H:246
barycentric.H
Foam::Vector< scalar >
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::tetrahedron::centre
Point centre() const
Return centre (centroid)
Definition: tetrahedronI.H:165
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:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::tetrahedron::Sd
vector Sd() const
Face area normal for side d.
Definition: tetrahedronI.H:158
Foam::tetrahedron::tetIntersectionList
FixedList< tetPoints, 200 > tetIntersectionList
Storage type for tets originating from intersecting tets.
Definition: tetrahedron.H:96
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Point
CGAL::Point_3< K > Point
Definition: CGALIndexedPolyhedron.H:53
rndGen
Random rndGen
Definition: createFields.H:23
Foam::tetrahedron::nEdges
Definition: tetrahedron.H:170
Foam::tetrahedron::storeOp::storeOp
storeOp(tetIntersectionList &, label &)
Definition: tetrahedronI.H:517
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
FixedList.H
Foam::tetrahedron::b
const Point & b() const
Definition: tetrahedronI.H:83
Foam::tetrahedron::bounds
boundBox bounds() const
Calculate the bounding box.
Definition: tetrahedron.C:342
tetrahedron.C
Foam::tetrahedron::Sc
vector Sc() const
Face area normal for side c.
Definition: tetrahedronI.H:151
Foam::tetrahedron::mag
scalar mag() const
Return volume.
Definition: tetrahedronI.H:172
tetrahedronI.H
Foam::tetrahedron::sliceWithPlane
void sliceWithPlane(const plane &pl, AboveTetOp &aboveOp, BelowTetOp &belowOp) const
Decompose tet into tets above and below plane.
Definition: tetrahedronI.H:1022
Foam::tetrahedron::barycentricToPoint
Point barycentricToPoint(const barycentric &bary) const
Calculate the point from the given barycentric coordinates.
Definition: tetrahedronI.H:256
Foam::tetrahedron::sumVolOp::operator()
void operator()(const tetPoints &)
Definition: tetrahedronI.H:507
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:65