tetCellI.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Description
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "IOstreams.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 :
36  FixedList<label, 4>(-1)
37 {}
38 
39 
41 (
42  const label a,
43  const label b,
44  const label c,
45  const label d
46 )
47 {
48  operator[](0) = a;
49  operator[](1) = b;
50  operator[](2) = c;
51  operator[](3) = d;
52 }
53 
54 
56 :
57  FixedList<label, 4>(lst)
58 {}
59 
60 
61 inline Foam::tetCell::tetCell(std::initializer_list<label> lst)
62 :
63  FixedList<label, 4>(lst)
64 {}
65 
66 
68 :
69  FixedList<label, 4>(is)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 inline Foam::triFace Foam::tetCell::face(const label facei) const
76 {
77  // Warning. Ordering of faces needs to be the same for a tetrahedron
78  // class, a tetrahedron cell shape model and a tetCell
79  static const label a[4] = {1, 0, 0, 0};
80  static const label b[4] = {2, 3, 1, 2};
81  static const label c[4] = {3, 2, 3, 1};
82 
83  #ifdef FULLDEBUG
84  if (facei < 0 || facei >= 4)
85  {
87  << "index out of range 0 -> 3. facei = " << facei
88  << abort(FatalError);
89  }
90  #endif
91 
92  return triFace
93  (
94  operator[](a[facei]),
95  operator[](b[facei]),
96  operator[](c[facei])
97  );
98 }
99 
100 
101 inline Foam::label Foam::tetCell::edgeFace(const label edgei) const
102 {
103  // Warning. Ordering of faces needs to be the same for a tetrahedron
104  // class, a tetrahedron cell shape model and a tetCell
105  static const label edgeFaces[6] = {2, 3, 1, 0, 0, 1};
106 
107  #ifdef FULLDEBUG
108  if (edgei < 0 || edgei >= 6)
109  {
111  << "edge index out of range 0 -> 5. edgei = " << edgei
112  << abort(FatalError);
113  }
114  #endif
115 
116  return edgeFaces[edgei];
117 }
118 
119 
121 (
122  const label edgei,
123  const label facei
124 ) const
125 {
126  // Warning. Ordering of faces needs to be the same for a tetrahedron
127  // class, a tetrahedron cell shape model and a tetCell
128  static const label adjacentFace[6][4] =
129  {
130  {-1, -1, 3, 2},
131  {-1, 3, -1, 1},
132  {-1, 2, 1, -1},
133  {2, -1, 0, -1},
134  {3, -1, -1, 0},
135  {1, 0, -1, -1}
136  };
137 
138  #ifdef FULLDEBUG
139  if (facei < 0 || facei >= 4)
140  {
142  << "face index out of range 0 -> 3. facei = " << facei
143  << abort(FatalError);
144  }
145 
146  if (edgei < 0 || edgei >= 6)
147  {
149  << "edge index out of range 0 -> 5. edgei = " << edgei
150  << abort(FatalError);
151  }
152  #endif
153 
154  return adjacentFace[edgei][facei];
155 }
156 
157 
158 inline Foam::edge Foam::tetCell::tetEdge(const label edgei) const
159 {
160  // Warning. Ordering of edges needs to be the same for a tetrahedron
161  // class, a tetrahedron cell shape model and a tetCell
162  static const label pt0[] = {0, 0, 0, 3, 1, 3};
163  static const label pt1[] = {1, 2, 3, 1, 2, 2};
164 
165  #ifdef FULLDEBUG
166  if (edgei < 0 || edgei >= 6)
167  {
169  << "index out of range 0 -> 5. edgei = " << edgei
170  << abort(FatalError);
171  }
172  #endif
173 
174  return edge(operator[](pt0[edgei]), operator[](pt1[edgei]));
175 }
176 
177 
179 {
180  return tetPointRef
181  (
182  points[operator[](0)],
183  points[operator[](1)],
184  points[operator[](2)],
185  points[operator[](3)]
186  );
187 }
188 
189 
190 // ************************************************************************* //
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::tetCell::tetCell
tetCell()
Construct null, with invalid point labels (-1)
Definition: tetCellI.H:34
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
Definition: tetPointRef.H:46
Foam::label
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:62
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::tetCell::edgeAdjacentFace
label edgeAdjacentFace(const label edgei, const label facei) const
Return face adjacent to the given face sharing the same edge.
Definition: tetCellI.H:121
Foam::tetCell::face
triFace face(const label facei) const
Return i-th face.
Definition: tetCellI.H:75
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::tetCell::tetEdge
edge tetEdge(const label edgei) const
Return i-th edge.
Definition: tetCellI.H:158
Foam::tetCell::edgeFace
label edgeFace(const label edgei) const
Return first face adjacent to the given edge.
Definition: tetCellI.H:101
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:70
Foam::tetCell::tet
tetPointRef tet(const UList< point > &points) const
Return the tetrahedron.
Definition: tetCellI.H:178
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::UList
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:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
triFace
face triFace(3)
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:65