Loading...
Searching...
No Matches
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-2022 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 Foam_tetrahedron_H
43#define Foam_tetrahedron_H
44
45#include "point.H"
46#include "pointHit.H"
47#include "primitiveFieldsFwd.H"
48#include "Random.H"
49#include "FixedList.H"
50#include "UList.H"
51#include "triangle.H"
52#include "treeBoundBox.H"
53#include "barycentric.H"
54
55// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57namespace Foam
58{
59
60// Forward Declarations
61class plane;
62
63template<class Point, class PointRef> class tetrahedron;
65template<class Point, class PointRef>
67
68template<class Point, class PointRef>
70
71
72// Common Typedefs
73
74//- A tetrahedron using referred points
76
77
78/*---------------------------------------------------------------------------*\
79 Class tetPoints Declaration
80\*---------------------------------------------------------------------------*/
81
82//- Tet point storage. Default constructable (tetrahedron is not)
83class tetPoints
84:
85 public FixedList<point, 4>
86{
87public:
88
89 // Generated Methods
90
91 //- Default construct
92 tetPoints() = default;
94 //- The front() accessor (from FixedList) has no purpose
95 void front() = delete;
96
97 //- The back() accessor (from FixedList) has no purpose
98 void back() = delete;
99
100
101 // Constructors
102
103 //- Construct from four points
104 inline tetPoints
105 (
106 const point& p0,
107 const point& p1,
108 const point& p2,
109 const point& p3
110 );
111
112 //- Construct from point references
113 inline explicit tetPoints(const tetPointRef& pts);
114
115 //- Construct from four points
116 inline tetPoints(const FixedList<point, 4>& pts);
117
118 //- Copy construct from subset of points
119 inline tetPoints
120 (
121 const UList<point>& points,
122 const FixedList<label, 4>& indices
123 );
124
125
126 // Member Functions
127
128 //- The first vertex
129 const point& a() const noexcept { return get<0>(); }
130
131 //- The second vertex
132 const point& b() const noexcept { return get<1>(); }
133
134 //- The third vertex
135 const point& c() const noexcept { return get<2>(); }
136
137 //- The fourth vertex
138 const point& d() const noexcept { return get<3>(); }
139
140 //- The first vertex
141 point& a() noexcept { return get<0>(); }
142
143 //- The second vertex
144 point& b() noexcept { return get<1>(); }
145
146 //- The third vertex
147 point& c() noexcept { return get<2>(); }
148
149 //- The fourth vertex
150 point& d() noexcept { return get<3>(); }
151
152 //- Return as tetrahedron reference
153 inline tetPointRef tet() const;
155 //- Invert tetrahedron by swapping third and fourth vertices
156 inline void flip();
157
158 //- The enclosing (bounding) box for the tetrahedron
159 inline Pair<point> box() const;
160
161 //- The bounding box for the tetrahedron
162 inline treeBoundBox bounds() const;
163};
165
166/*---------------------------------------------------------------------------*\
167 Class tetrahedron Declaration
168\*---------------------------------------------------------------------------*/
170template<class Point, class PointRef>
171class tetrahedron
172{
173public:
175 // Public Typedefs
176
177 //- The point type
178 typedef Point point_type;
180 //- Storage type for tets originating from intersecting tets.
181 // (can possibly be smaller than 200)
183
184
185 // Classes for use in sliceWithPlane.
186 // What to do with decomposition of tetrahedron.
187
188 //- Dummy
189 struct dummyOp
190 {
191 void operator()(const tetPoints&) const noexcept
192 {}
193 };
194
195 //- Sum resulting volumes
196 struct sumVolOp
197 {
198 scalar vol_;
199
200 constexpr sumVolOp() noexcept : vol_(0) {}
201
202 scalar total() const noexcept { return vol_; }
203
204 inline void operator()(const tetPoints&);
205 };
206
207 //- Store resulting tets
208 struct storeOp
209 {
211 label& nTets_;
212
213 inline storeOp(tetIntersectionList&, label&);
214
215 inline void operator()(const tetPoints&);
216 };
218private:
219
220 // Private Data
221
222 PointRef a_, b_, c_, d_;
223
225 // Private Member Functions
226
227 inline static point planeIntersection
228 (
230 const tetPoints&,
231 const label,
232 const label
233 );
234
235 template<class TetOp>
236 inline static void decomposePrism
237 (
239 TetOp& op
240 );
241
242 template<class AboveTetOp, class BelowTetOp>
243 inline static void tetSliceWithPlane
245 const plane& pl,
246 const tetPoints& tet,
247 AboveTetOp& aboveOp,
248 BelowTetOp& belowOp
249 );
250
251
252public:
253
254 // Member Constants
255
256 enum
257 {
258 nVertices = 4, // Number of vertices in tetrahedron
259 nEdges = 6 // Number of edges in tetrahedron
260 };
261
262
263 // Constructors
264
265 //- Construct from four points
266 inline tetrahedron
267 (
268 const Point& p0,
269 const Point& p1,
270 const Point& p2,
271 const Point& p3
272 );
273
274 //- Construct from four points
275 inline tetrahedron(const FixedList<Point, 4>& pts);
276
277 //- Construct from four points in the list of points
278 inline tetrahedron
279 (
280 const UList<Point>& points,
281 const FixedList<label, 4>& indices
282 );
283
284 //- Construct from Istream
285 inline explicit tetrahedron(Istream&);
286
287
288 // Member Functions
289
290 // Access
291
292 //- Return vertex a
293 const Point& a() const noexcept { return a_; }
294
295 //- Return vertex b
296 const Point& b() const noexcept { return b_; }
297
298 //- Return vertex c
299 const Point& c() const noexcept { return c_; }
300
301 //- Return vertex d
302 const Point& d() const noexcept { return d_; }
303
304 //- Return i-th face
305 inline triPointRef tri(const label facei) const;
308 // Tet properties (static calculations)
309
310 //- The enclosing (bounding) box for four points
311 inline static Pair<Point> box
312 (
313 const Point& p0,
314 const Point& p1,
315 const Point& p2,
316 const Point& p3
317 );
318
319
320 // Properties
321
322 //- Face area normal for side a
323 inline vector Sa() const;
324
325 //- Face area normal for side b
326 inline vector Sb() const;
328 //- Face area normal for side c
329 inline vector Sc() const;
330
331 //- Face area normal for side d
332 inline vector Sd() const;
333
334 //- Return centre (centroid)
335 inline Point centre() const;
336
337 //- Return volume
338 inline scalar mag() const;
339
340 //- The enclosing (bounding) box for the tetrahedron
341 inline Pair<Point> box() const;
342
343 //- The bounding box for the tetrahedron
344 inline treeBoundBox bounds() const;
345
346 //- Return circum-centre
347 inline Point circumCentre() const;
348
349 //- Return circum-radius
350 inline scalar circumRadius() const;
352 //- Return quality: Ratio of tetrahedron and circum-sphere
353 //- volume, scaled so that a regular tetrahedron has a
354 //- quality of 1
355 inline scalar quality() const;
357 //- Return a random point in the tetrahedron from a
358 //- uniform distribution
359 inline Point randomPoint(Random& rndGen) const;
360
361 //- Calculate the point from the given barycentric coordinates.
362 inline Point barycentricToPoint(const barycentric& bary) const;
363
364 //- Calculate the barycentric coordinates from the given point
365 inline barycentric pointToBarycentric(const point& pt) const;
367 //- Calculate the barycentric coordinates from the given point.
368 // Returns the determinant.
369 inline scalar pointToBarycentric
370 (
371 const point& pt,
372 barycentric& bary
373 ) const;
374
375 //- Return nearest point to p on tetrahedron. Is p itself
376 // if inside.
377 inline pointHit nearestPoint(const point& p) const;
378
379 //- Return true if point is inside tetrahedron
380 inline bool inside(const point& pt) const;
381
382 //- Decompose tet into tets above and below plane
383 template<class AboveTetOp, class BelowTetOp>
384 inline void sliceWithPlane
385 (
386 const plane& pl,
387 AboveTetOp& aboveOp,
388 BelowTetOp& belowOp
389 ) const;
390
391 //- Decompose tet into tets inside and outside other tet
392 inline void tetOverlap
395 tetIntersectionList& insideTets,
396 label& nInside,
397 tetIntersectionList& outsideTets,
398 label& nOutside
399 ) const;
400
401
402 //- Return (min)containment sphere, i.e. the smallest sphere with
403 // all points inside. Returns pointHit with:
404 // - hit : if sphere is equal to circumsphere
405 // (biggest sphere)
406 // - point : centre of sphere
407 // - distance : radius of sphere
408 // - eligiblemiss: false
409 // Tol (small compared to 1, e.g. 1e-9) is used to determine
410 // whether point is inside: mag(pt - ctr) < (1+tol)*radius.
411 pointHit containmentSphere(const scalar tol) const;
412
413 //- Fill buffer with shape function products
414 void gradNiSquared(scalarField& buffer) const;
415
416 void gradNiDotGradNj(scalarField& buffer) const;
417
418 void gradNiGradNi(tensorField& buffer) const;
419
420 void gradNiGradNj(tensorField& buffer) const;
421
422
423 // IOstream Operators
424
425 friend Istream& operator>> <Point, PointRef>
426 (
427 Istream&,
429 );
430
432 (
434 const tetrahedron&
435 );
436};
437
439// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
440
441} // End namespace Foam
442
443// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
444
446
447// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448
449#ifdef NoRepository
450 #include "tetrahedron.C"
451#endif
452
453// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
454
455#endif
457// ************************************************************************* //
CGAL::Point_3< K > Point
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
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
Tet point storage. Default constructable (tetrahedron is not).
Definition tetrahedron.H:85
void flip()
Invert tetrahedron by swapping third and fourth vertices.
point & b() noexcept
The second vertex.
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.
point & a() noexcept
The first vertex.
const point & d() const noexcept
The fourth vertex.
const point & b() const noexcept
The second vertex.
Pair< point > box() const
The enclosing (bounding) box for the tetrahedron.
point & d() noexcept
The fourth vertex.
tetPointRef tet() const
Return as tetrahedron reference.
point & c() noexcept
The third vertex.
const point & a() const noexcept
The first vertex.
tetPoints()=default
Default construct.
treeBoundBox bounds() const
The bounding box for the tetrahedron.
A tetrahedron primitive.
Point circumCentre() const
Return circum-centre.
barycentric pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
pointHit containmentSphere(const scalar tol) const
Return (min)containment sphere, i.e. the smallest sphere with.
Definition tetrahedron.C:28
scalar pointToBarycentric(const point &pt, barycentric &bary) const
Calculate the barycentric coordinates from the given point.
void gradNiDotGradNj(scalarField &buffer) const
void gradNiGradNi(tensorField &buffer) const
vector Sd() const
Face area normal for side d.
triPointRef tri(const label facei) const
Return i-th face.
tetrahedron(const FixedList< Point, 4 > &pts)
Construct from four points.
scalar circumRadius() const
Return circum-radius.
Point centre() const
Return centre (centroid).
vector Sc() const
Face area normal for side c.
bool inside(const point &pt) const
Return true if point is inside tetrahedron.
vector Sb() const
Face area normal for side b.
tetrahedron(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Construct from four points.
const Point & b() const noexcept
Return vertex b.
scalar mag() const
Return volume.
friend Ostream & operator(Ostream &, const tetrahedron &)
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
pointHit nearestPoint(const point &p) const
Return nearest point to p on tetrahedron. Is p itself.
Point randomPoint(Random &rndGen) const
Return a random point in the tetrahedron from a uniform distribution.
tetrahedron(const UList< Point > &points, const FixedList< label, 4 > &indices)
Construct from four points in the list of points.
Pair< Point > box() const
The enclosing (bounding) box for the tetrahedron.
FixedList< tetPoints, 200 > tetIntersectionList
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.
treeBoundBox bounds() const
The bounding box for the tetrahedron.
tetrahedron(Istream &)
Construct from Istream.
scalar quality() const
Return quality: Ratio of tetrahedron and circum-sphere volume, scaled so that a regular tetrahedron h...
const Point & c() const noexcept
Return vertex c.
void gradNiSquared(scalarField &buffer) const
Fill buffer with shape function products.
const Point & a() const noexcept
Return vertex a.
vector Sa() const
Face area normal for side a.
const Point & d() const noexcept
Return vertex d.
static Pair< Point > box(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
The enclosing (bounding) box for four points.
Standard boundBox with extra functionality for use in octree.
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
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition tetrahedron.H:72
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Istream & operator>>(Istream &, directionInfo &)
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
const pointField & pts
void operator()(const tetPoints &) const noexcept
Store resulting tets.
tetIntersectionList & tets_
void operator()(const tetPoints &)
storeOp(tetIntersectionList &, label &)
void operator()(const tetPoints &)
constexpr sumVolOp() noexcept
scalar total() const noexcept
Random rndGen