00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include <OpenFOAM/triangle.H>
00027 #include <OpenFOAM/IOstreams.H>
00028
00029
00030
00031 namespace Foam
00032 {
00033
00034
00035
00036 template<class Point, class PointRef>
00037 inline tetrahedron<Point, PointRef>::tetrahedron
00038 (
00039 const Point& a,
00040 const Point& b,
00041 const Point& c,
00042 const Point& d
00043 )
00044 :
00045 a_(a),
00046 b_(b),
00047 c_(c),
00048 d_(d)
00049 {}
00050
00051
00052 template<class Point, class PointRef>
00053 inline tetrahedron<Point, PointRef>::tetrahedron(Istream& is)
00054 {
00055
00056 is.readBegin("tetrahedron");
00057
00058 is >> a_ >> b_ >> c_ >> d_;
00059
00060
00061 is.readEnd("tetrahedron");
00062
00063
00064 is.check("tetrahedron::tetrahedron(Istream& is)");
00065 }
00066
00067
00068
00069
00070 template<class Point, class PointRef>
00071 inline const Point& tetrahedron<Point, PointRef>::a() const
00072 {
00073 return a_;
00074 }
00075
00076
00077 template<class Point, class PointRef>
00078 inline const Point& tetrahedron<Point, PointRef>::b() const
00079 {
00080 return b_;
00081 }
00082
00083
00084 template<class Point, class PointRef>
00085 inline const Point& tetrahedron<Point, PointRef>::c() const
00086 {
00087 return c_;
00088 }
00089
00090
00091 template<class Point, class PointRef>
00092 inline const Point& tetrahedron<Point, PointRef>::d() const
00093 {
00094 return d_;
00095 }
00096
00097
00098 template<class Point, class PointRef>
00099 inline vector tetrahedron<Point, PointRef>::Sa() const
00100 {
00101 return triangle<Point, PointRef>(b_, c_, d_).normal();
00102 }
00103
00104
00105 template<class Point, class PointRef>
00106 inline vector tetrahedron<Point, PointRef>::Sb() const
00107 {
00108 return triangle<Point, PointRef>(a_, d_, c_).normal();
00109 }
00110
00111
00112 template<class Point, class PointRef>
00113 inline vector tetrahedron<Point, PointRef>::Sc() const
00114 {
00115 return triangle<Point, PointRef>(a_, b_, d_).normal();
00116 }
00117
00118
00119 template<class Point, class PointRef>
00120 inline vector tetrahedron<Point, PointRef>::Sd() const
00121 {
00122 return triangle<Point, PointRef>(a_, c_, b_).normal();
00123 }
00124
00125
00126 template<class Point, class PointRef>
00127 inline Point tetrahedron<Point, PointRef>::centre() const
00128 {
00129 return 0.25*(a_ + b_ + c_ + d_);
00130 }
00131
00132
00133 template<class Point, class PointRef>
00134 inline scalar tetrahedron<Point, PointRef>::mag() const
00135 {
00136 return (1.0/6.0)*(((b_ - a_) ^ (c_ - a_)) & (d_ - a_));
00137 }
00138
00139
00140 template<class Point, class PointRef>
00141 inline Point tetrahedron<Point, PointRef>::circumCentre() const
00142 {
00143 vector a = b_ - a_;
00144 vector b = c_ - a_;
00145 vector c = d_ - a_;
00146
00147 scalar lamda = magSqr(c) - (a&c);
00148 scalar mu = magSqr(b) - (a&b);
00149
00150 vector ba = b^a;
00151 vector ca = c^a;
00152
00153 return a_ + 0.5*(a + (lamda*ba - mu*ca)/(c&ba));
00154 }
00155
00156
00157 template<class Point, class PointRef>
00158 inline scalar tetrahedron<Point, PointRef>::circumRadius() const
00159 {
00160 return Foam::mag(a_ - circumCentre());
00161 }
00162
00163
00164
00165
00166 template<class point, class pointRef>
00167 inline Istream& operator>>(Istream& is, tetrahedron<point, pointRef>& t)
00168 {
00169
00170 is.readBegin("tetrahedron");
00171
00172 is >> t.a_ >> t.b_ >> t.c_ >> t.d_;
00173
00174
00175 is.readEnd("tetrahedron");
00176
00177
00178 is.check("Istream& operator>>(Istream&, tetrahedron&)");
00179
00180 return is;
00181 }
00182
00183
00184 template<class Point, class PointRef>
00185 inline Ostream& operator<<(Ostream& os, const tetrahedron<Point, PointRef>& t)
00186 {
00187 os << nl
00188 << token::BEGIN_LIST
00189 << t.a_ << token::SPACE << t.b_
00190 << token::SPACE << t.c_ << token::SPACE << t.d_
00191 << token::END_LIST;
00192
00193 return os;
00194 }
00195
00196
00197
00198
00199 }
00200
00201