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 -------------------------------------------------------------------------------
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 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 
57 namespace Foam
58 {
59 
60 // Forward Declarations
61 class plane;
62 
63 template<class Point, class PointRef> class tetrahedron;
64 
65 template<class Point, class PointRef>
67 
68 template<class Point, class PointRef>
69 inline Ostream& operator<<(Ostream&, const tetrahedron<Point, 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)
83 class tetPoints
84 :
85  public FixedList<point, 4>
86 {
87 public:
88 
89  // Generated Methods
90 
91  //- Default construct
92  tetPoints() = default;
93 
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 \*---------------------------------------------------------------------------*/
170 template<class Point, class PointRef>
171 class tetrahedron
172 {
173 public:
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. What to do with decomposition
186  // of tet.
187 
188  //- Dummy
189  class dummyOp
190  {
191  public:
192  inline void operator()(const tetPoints&);
193  };
194 
195  //- Sum resulting volumes
196  class sumVolOp
197  {
198  public:
199  scalar vol_;
200 
201  inline sumVolOp();
202 
203  inline void operator()(const tetPoints&);
204  };
205 
206  //- Store resulting tets
207  class storeOp
208  {
209  tetIntersectionList& tets_;
210  label& nTets_;
211 
212  public:
213  inline storeOp(tetIntersectionList&, label&);
214 
215  inline void operator()(const tetPoints&);
216  };
218 private:
219 
220  // Private Data
221 
222  PointRef a_, b_, c_, d_;
223 
225  // Private Member Functions
226 
227  inline static point planeIntersection
228  (
229  const FixedList<scalar, 4>&,
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
244  (
245  const plane& pl,
246  const tetPoints& tet,
247  AboveTetOp& aboveOp,
248  BelowTetOp& belowOp
249  );
250 
251 
252 public:
253 
254  // Member Constants
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;
327 
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
393  (
394  const tetrahedron<Point, PointRef>& tetB,
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&,
428  tetrahedron&
429  );
430 
431  friend Ostream& operator<< <Point, PointRef>
432  (
433  Ostream&,
434  const tetrahedron&
435  );
436 };
437 
438 
439 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
440 
441 } // End namespace Foam
442 
443 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
444 
445 #include "tetrahedronI.H"
446 
447 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448 
449 #ifdef NoRepository
450  #include "tetrahedron.C"
451 #endif
452 
453 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
454 
455 #endif
456 
457 // ************************************************************************* //
treeBoundBox bounds() const
The bounding box for the tetrahedron.
Definition: tetrahedronI.H:174
A tetrahedron primitive.
Definition: tetrahedron.H:58
const Point & c() const noexcept
Return vertex c.
Definition: tetrahedron.H:361
void operator()(const tetPoints &)
Definition: tetrahedronI.H:593
Pair< Point > box() const
The enclosing (bounding) box for the tetrahedron.
Definition: tetrahedronI.H:155
pointHit containmentSphere(const scalar tol) const
Return (min)containment sphere, i.e. the smallest sphere with.
Definition: tetrahedron.C:28
void operator()(const tetPoints &)
Definition: tetrahedronI.H:615
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
scalar circumRadius() const
Return circum-radius.
Definition: tetrahedronI.H:286
Point circumCentre() const
Return circum-centre.
Definition: tetrahedronI.H:259
tetrahedron(const Point &p0, const Point &p1, const Point &p2, const Point &p3)
Construct from four points.
Definition: tetrahedronI.H:70
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetrahedron.H:72
tetPointRef tet() const
Return as tetrahedron reference.
Definition: tetrahedronI.H:120
vector Sd() const
Face area normal for side d.
Definition: tetrahedronI.H:238
const point & a() const noexcept
The first vertex.
Definition: tetrahedron.H:144
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void gradNiSquared(scalarField &buffer) const
Fill buffer with shape function products.
Definition: tetrahedron.C:240
Random rndGen
Definition: createFields.H:23
void operator()(const tetPoints &)
Definition: tetrahedronI.H:578
storeOp(tetIntersectionList &, label &)
Definition: tetrahedronI.H:603
tetPoints()=default
Default construct.
const point & c() const noexcept
The third vertex.
Definition: tetrahedron.H:154
FixedList< tetPoints, 200 > tetIntersectionList
Storage type for tets originating from intersecting tets.
Definition: tetrahedron.H:224
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:90
Point barycentricToPoint(const barycentric &bary) const
Calculate the point from the given barycentric coordinates.
Definition: tetrahedronI.H:336
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
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.
const pointField & points
const Point & a() const noexcept
Return vertex a.
Definition: tetrahedron.H:351
scalar quality() const
Return quality: Ratio of tetrahedron and circum-sphere volume, scaled so that a regular tetrahedron h...
Definition: tetrahedronI.H:312
Istream & operator>>(Istream &, directionInfo &)
void gradNiGradNi(tensorField &buffer) const
Definition: tetrahedron.C:291
vector Sb() const
Face area normal for side b.
Definition: tetrahedronI.H:224
scalar mag() const
Return volume.
Definition: tetrahedronI.H:252
vector Sa() const
Face area normal for side a.
Definition: tetrahedronI.H:217
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
Point centre() const
Return centre (centroid)
Definition: tetrahedronI.H:245
bool inside(const point &pt) const
Return true if point is inside tetrahedron.
Definition: tetrahedronI.H:515
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:105
Random number generator.
Definition: Random.H:55
Pair< point > box() const
The enclosing (bounding) box for the tetrahedron.
Definition: tetrahedronI.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
const Point & b() const noexcept
Return vertex b.
Definition: tetrahedron.H:356
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
void back()=delete
The back() accessor (from FixedList) has no purpose.
void sliceWithPlane(const plane &pl, AboveTetOp &aboveOp, BelowTetOp &belowOp) const
Decompose tet into tets above and below plane.
CGAL::Point_3< K > Point
const point & b() const noexcept
The second vertex.
Definition: tetrahedron.H:149
triPointRef tri(const label facei) const
Return i-th face.
Definition: tetrahedronI.H:185
barycentric pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: tetrahedronI.H:346
Tet point storage. Default constructable (tetrahedron is not)
Definition: tetrahedron.H:82
void gradNiGradNj(tensorField &buffer) const
Definition: tetrahedron.C:309
pointHit nearestPoint(const point &p) const
Return nearest point to p on tetrahedron. Is p itself.
Definition: tetrahedronI.H:401
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
const Point & d() const noexcept
Return vertex d.
Definition: tetrahedron.H:366
treeBoundBox bounds() const
The bounding box for the tetrahedron.
Definition: tetrahedronI.H:168
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
Point point_type
The point type.
Definition: tetrahedron.H:217
void gradNiDotGradNj(scalarField &buffer) const
Definition: tetrahedron.C:261
Point randomPoint(Random &rndGen) const
Return a random point in the tetrahedron from a uniform distribution.
Definition: tetrahedronI.H:326
volScalarField & p
vector Sc() const
Face area normal for side c.
Definition: tetrahedronI.H:231
const point & d() const noexcept
The fourth vertex.
Definition: tetrahedron.H:159
void flip()
Invert tetrahedron by swapping third and fourth vertices.
Definition: tetrahedronI.H:126
const volScalarField & p0
Definition: EEqn.H:36
Namespace for OpenFOAM.
const pointField & pts
void front()=delete
The front() accessor (from FixedList) has no purpose.