initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / OpenFOAM / meshes / primitiveShapes / tetrahedron / tetrahedronI.H
blobc6b08b54ab1f0d865ce6245a24bad83e3c2bda63
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Description
27 \*---------------------------------------------------------------------------*/
29 #include "triangle.H"
30 #include "IOstreams.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
39 template<class Point, class PointRef>
40 inline tetrahedron<Point, PointRef>::tetrahedron
42     const Point& a,
43     const Point& b,
44     const Point& c,
45     const Point& d
48     a_(a),
49     b_(b),
50     c_(c),
51     d_(d)
55 template<class Point, class PointRef>
56 inline tetrahedron<Point, PointRef>::tetrahedron(Istream& is)
58     // Read beginning of tetrahedron point pair
59     is.readBegin("tetrahedron");
61     is >> a_ >> b_ >> c_ >> d_;
63     // Read end of tetrahedron point pair
64     is.readEnd("tetrahedron");
66     // Check state of Istream
67     is.check("tetrahedron::tetrahedron(Istream& is)");
71 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
73 template<class Point, class PointRef>
74 inline const Point& tetrahedron<Point, PointRef>::a() const
76     return a_;
80 template<class Point, class PointRef>
81 inline const Point& tetrahedron<Point, PointRef>::b() const
83     return b_;
87 template<class Point, class PointRef>
88 inline const Point& tetrahedron<Point, PointRef>::c() const
90     return c_;
94 template<class Point, class PointRef>
95 inline const Point& tetrahedron<Point, PointRef>::d() const
97     return d_;
101 template<class Point, class PointRef>
102 inline vector tetrahedron<Point, PointRef>::Sa() const
104     return triangle<Point, PointRef>(b_, c_, d_).normal();
108 template<class Point, class PointRef>
109 inline vector tetrahedron<Point, PointRef>::Sb() const
111     return triangle<Point, PointRef>(a_, d_, c_).normal();
115 template<class Point, class PointRef>
116 inline vector tetrahedron<Point, PointRef>::Sc() const
118     return triangle<Point, PointRef>(a_, b_, d_).normal();
122 template<class Point, class PointRef>
123 inline vector tetrahedron<Point, PointRef>::Sd() const
125     return triangle<Point, PointRef>(a_, c_, b_).normal();
129 template<class Point, class PointRef>
130 inline scalar tetrahedron<Point, PointRef>::mag() const
132     return (1.0/6.0)*(((b_ - a_) ^ (c_ - a_)) & (d_ - a_));
136 template<class Point, class PointRef>
137 inline vector tetrahedron<Point, PointRef>::circumCentre() const
139     vector a = b_ - a_;
140     vector b = c_ - a_;
141     vector c = d_ - a_;
143     scalar lamda = magSqr(c) - (a&c);
144     scalar mu = magSqr(b) - (a&b);
146     vector ba = b^a;
147     vector ca = c^a;
149     return a_ + 0.5*(a + (lamda*ba - mu*ca)/(c&ba));
153 template<class Point, class PointRef>
154 inline scalar tetrahedron<Point, PointRef>::circumRadius() const
156     return Foam::mag(a_ - circumCentre());
160 // * * * * * * * * * * * * * * * Ostream Operator  * * * * * * * * * * * * * //
162 template<class point, class pointRef>
163 inline Istream& operator>>(Istream& is, tetrahedron<point, pointRef>& t)
165     // Read beginning of tetrahedron point pair
166     is.readBegin("tetrahedron");
168     is >> t.a_ >> t.b_ >> t.c_ >> t.d_;
170     // Read end of tetrahedron point pair
171     is.readEnd("tetrahedron");
173     // Check state of Ostream
174     is.check("Istream& operator>>(Istream&, tetrahedron&)");
176     return is;
180 template<class Point, class PointRef>
181 inline Ostream& operator<<(Ostream& os, const tetrahedron<Point, PointRef>& t)
183     os  << nl
184         << token::BEGIN_LIST
185         << t.a_ << token::SPACE << t.b_
186         << token::SPACE << t.c_ << token::SPACE << t.d_
187         << token::END_LIST;
189     return os;
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 } // End namespace Foam
197 // ************************************************************************* //