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