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