FreeFOAM The Cross-Platform CFD Toolkit
tetrahedronI.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) 1991-2010 OpenCFD Ltd.
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 \*---------------------------------------------------------------------------*/
25 
26 #include <OpenFOAM/triangle.H>
27 #include <OpenFOAM/IOstreams.H>
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Point, class PointRef>
38 (
39  const Point& a,
40  const Point& b,
41  const Point& c,
42  const Point& d
43 )
44 :
45  a_(a),
46  b_(b),
47  c_(c),
48  d_(d)
49 {}
50 
51 
52 template<class Point, class PointRef>
54 {
55  // Read beginning of tetrahedron point pair
56  is.readBegin("tetrahedron");
57 
58  is >> a_ >> b_ >> c_ >> d_;
59 
60  // Read end of tetrahedron point pair
61  is.readEnd("tetrahedron");
62 
63  // Check state of Istream
64  is.check("tetrahedron::tetrahedron(Istream& is)");
65 }
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class Point, class PointRef>
71 inline const Point& tetrahedron<Point, PointRef>::a() const
72 {
73  return a_;
74 }
75 
76 
77 template<class Point, class PointRef>
78 inline const Point& tetrahedron<Point, PointRef>::b() const
79 {
80  return b_;
81 }
82 
83 
84 template<class Point, class PointRef>
85 inline const Point& tetrahedron<Point, PointRef>::c() const
86 {
87  return c_;
88 }
89 
90 
91 template<class Point, class PointRef>
92 inline const Point& tetrahedron<Point, PointRef>::d() const
93 {
94  return d_;
95 }
96 
97 
98 template<class Point, class PointRef>
100 {
101  return triangle<Point, PointRef>(b_, c_, d_).normal();
102 }
103 
104 
105 template<class Point, class PointRef>
107 {
108  return triangle<Point, PointRef>(a_, d_, c_).normal();
109 }
110 
111 
112 template<class Point, class PointRef>
114 {
115  return triangle<Point, PointRef>(a_, b_, d_).normal();
116 }
117 
118 
119 template<class Point, class PointRef>
121 {
122  return triangle<Point, PointRef>(a_, c_, b_).normal();
123 }
124 
125 
126 template<class Point, class PointRef>
128 {
129  return 0.25*(a_ + b_ + c_ + d_);
130 }
131 
132 
133 template<class Point, class PointRef>
135 {
136  return (1.0/6.0)*(((b_ - a_) ^ (c_ - a_)) & (d_ - a_));
137 }
138 
139 
140 template<class Point, class PointRef>
142 {
143  vector a = b_ - a_;
144  vector b = c_ - a_;
145  vector c = d_ - a_;
146 
147  scalar lamda = magSqr(c) - (a&c);
148  scalar mu = magSqr(b) - (a&b);
149 
150  vector ba = b^a;
151  vector ca = c^a;
152 
153  return a_ + 0.5*(a + (lamda*ba - mu*ca)/(c&ba));
154 }
155 
156 
157 template<class Point, class PointRef>
159 {
160  return Foam::mag(a_ - circumCentre());
161 }
162 
163 
164 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
165 
166 template<class point, class pointRef>
168 {
169  // Read beginning of tetrahedron point pair
170  is.readBegin("tetrahedron");
171 
172  is >> t.a_ >> t.b_ >> t.c_ >> t.d_;
173 
174  // Read end of tetrahedron point pair
175  is.readEnd("tetrahedron");
176 
177  // Check state of Ostream
178  is.check("Istream& operator>>(Istream&, tetrahedron&)");
179 
180  return is;
181 }
182 
183 
184 template<class Point, class PointRef>
185 inline Ostream& operator<<(Ostream& os, const tetrahedron<Point, PointRef>& t)
186 {
187  os << nl
189  << t.a_ << token::SPACE << t.b_
190  << token::SPACE << t.c_ << token::SPACE << t.d_
191  << token::END_LIST;
192 
193  return os;
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace Foam
200 
201 // ************************ vim: set sw=4 sts=4 et: ************************ //