Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / OpenFOAM / meshes / primitiveShapes / tetrahedron / tetrahedron.H
blob375da8525a58fb62d9c4dce7be3c6a7cd5618e73
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2011 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Class
25     Foam::tetrahedron
27 Description
28     A tetrahedron primitive.
30     Ordering of edges needs to be the same for a tetrahedron
31     class, a tetrahedron cell shape model and a tetCell.
33 SourceFiles
34     tetrahedronI.H
35     tetrahedron.C
37 \*---------------------------------------------------------------------------*/
39 #ifndef tetrahedron_H
40 #define tetrahedron_H
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"
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 namespace Foam
55 class Istream;
56 class Ostream;
58 // Forward declaration of friend functions and operators
60 template<class Point, class PointRef> class tetrahedron;
62 template<class Point, class PointRef>
63 inline Istream& operator>>
65     Istream&,
66     tetrahedron<Point, PointRef>&
69 template<class Point, class PointRef>
70 inline Ostream& operator<<
72     Ostream&,
73     const tetrahedron<Point, PointRef>&
77 /*---------------------------------------------------------------------------*\
78                            class tetrahedron Declaration
79 \*---------------------------------------------------------------------------*/
81 template<class Point, class PointRef>
82 class tetrahedron
84     // Private data
86         PointRef a_, b_, c_, d_;
89 public:
91     // Member constants
93         enum
94         {
95             nVertices = 4,  // Number of vertices in tetrahedron
96             nEdges = 6      // Number of edges in tetrahedron
97         };
100     // Constructors
102         //- Construct from points
103         inline tetrahedron
104         (
105             const Point& a,
106             const Point& b,
107             const Point& c,
108             const Point& d
109         );
111         //- Construct from four points in the list of points
112         inline tetrahedron
113         (
114             const UList<Point>&,
115             const FixedList<label, 4>& indices
116         );
118         //- Construct from Istream
119         inline tetrahedron(Istream&);
122     // Member Functions
124         // Access
126             //- Return vertices
127             inline const Point& a() const;
129             inline const Point& b() const;
131             inline const Point& c() const;
133             inline const Point& d() const;
136         // Properties
138             //- Return face normal
139             inline vector Sa() const;
141             inline vector Sb() const;
143             inline vector Sc() const;
145             inline vector Sd() const;
147             //- Return centre (centroid)
148             inline Point centre() const;
150             //- Return volume
151             inline scalar mag() const;
153             //- Return circum-centre
154             inline Point circumCentre() const;
156             //- Return circum-radius
157             inline scalar circumRadius() const;
159             //- Return quality: Ratio of tetrahedron and circum-sphere
160             //  volume, scaled so that a regular tetrahedron has a
161             //  quality of 1
162             inline scalar quality() const;
164             //- Return a random point in the tetrahedron from a
165             //  uniform distribution
166             inline Point randomPoint(Random& rndGen) const;
168             //- Return a random point in the tetrahedron from a
169             //  uniform distribution
170             inline Point randomPoint(cachedRandom& rndGen) const;
172             //- Calculate the barycentric coordinates of the given
173             //  point, in the same order as a, b, c, d.  Returns the
174             //  determinant of the solution.
175             inline scalar barycentric
176             (
177                 const point& pt,
178                 List<scalar>& bary
179             ) const;
181             //- Return nearest point to p on tetrahedron
182             inline pointHit nearestPoint(const point& p) const;
184             //- Return true if point is inside tetrahedron
185             inline bool inside(const point& pt) const;
187             //- Return (min)containment sphere, i.e. the smallest sphere with
188             //  all points inside. Returns pointHit with:
189             //  - hit         : if sphere is equal to circumsphere
190             //                  (biggest sphere)
191             //  - point       : centre of sphere
192             //  - distance    : radius of sphere
193             //  - eligiblemiss: false
194             // Tol (small compared to 1, e.g. 1E-9) is used to determine
195             // whether point is inside: mag(pt - ctr) < (1+tol)*radius.
196             pointHit containmentSphere(const scalar tol) const;
198             //- Fill buffer with shape function products
199             void gradNiSquared(scalarField& buffer) const;
201             void gradNiDotGradNj(scalarField& buffer) const;
203             void gradNiGradNi(tensorField& buffer) const;
205             void gradNiGradNj(tensorField& buffer) const;
208     // IOstream operators
210         friend Istream& operator>> <Point, PointRef>
211         (
212             Istream&,
213             tetrahedron&
214         );
216         friend Ostream& operator<< <Point, PointRef>
217         (
218             Ostream&,
219             const tetrahedron&
220         );
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 } // End namespace Foam
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 #include "tetrahedronI.H"
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 #ifdef NoRepository
235 #   include "tetrahedron.C"
236 #endif
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 #endif
242 // ************************************************************************* //