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  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 :
33  FixedList<label, 4>(-1)
34 {}
35 
36 
38 (
39  const label a,
40  const label b,
41  const label c,
42  const label d
43 )
44 {
45  operator[](0) = a;
46  operator[](1) = b;
47  operator[](2) = c;
48  operator[](3) = d;
49 }
50 
51 
52 inline Foam::tetCell::tetCell(std::initializer_list<label> list)
53 :
54  FixedList<label, 4>(list)
55 {}
56 
57 
59 :
60  FixedList<label, 4>(list)
61 {}
62 
63 
65 (
66  const labelUList& list,
67  const FixedList<label, 4>& indices
68 )
69 :
70  FixedList<label, 4>(list, indices)
71 {}
72 
73 
75 :
76  FixedList<label, 4>(is)
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
82 inline Foam::triFace Foam::tetCell::face(const label facei) const
83 {
84  // Warning. Ordering of faces needs to be the same for a tetrahedron
85  // class, a tetrahedron cell shape model and a tetCell
86  static const label a[4] = {1, 0, 0, 0};
87  static const label b[4] = {2, 3, 1, 2};
88  static const label c[4] = {3, 2, 3, 1};
89 
90  #ifdef FULLDEBUG
91  if (facei < 0 || facei >= 4)
92  {
94  << "index out of range 0 -> 3. facei = " << facei
95  << abort(FatalError);
96  }
97  #endif
98 
99  return triFace
100  (
101  operator[](a[facei]),
102  operator[](b[facei]),
103  operator[](c[facei])
104  );
105 }
106 
107 
108 inline Foam::label Foam::tetCell::edgeFace(const label edgei) const
109 {
110  // Warning. Ordering of faces needs to be the same for a tetrahedron
111  // class, a tetrahedron cell shape model and a tetCell
112  static const label edgeFaces[6] = {2, 3, 1, 0, 0, 1};
113 
114  #ifdef FULLDEBUG
115  if (edgei < 0 || edgei >= 6)
116  {
118  << "edge index out of range 0 -> 5. edgei = " << edgei
119  << abort(FatalError);
120  }
121  #endif
122 
123  return edgeFaces[edgei];
124 }
125 
126 
127 inline Foam::label Foam::tetCell::edgeAdjacentFace
128 (
129  const label edgei,
130  const label facei
131 ) const
132 {
133  // Warning. Ordering of faces needs to be the same for a tetrahedron
134  // class, a tetrahedron cell shape model and a tetCell
135  static const label adjacentFace[6][4] =
136  {
137  {-1, -1, 3, 2},
138  {-1, 3, -1, 1},
139  {-1, 2, 1, -1},
140  {2, -1, 0, -1},
141  {3, -1, -1, 0},
142  {1, 0, -1, -1}
143  };
144 
145  #ifdef FULLDEBUG
146  if (facei < 0 || facei >= 4)
147  {
149  << "face index out of range 0 -> 3. facei = " << facei
150  << abort(FatalError);
151  }
152 
153  if (edgei < 0 || edgei >= 6)
154  {
156  << "edge index out of range 0 -> 5. edgei = " << edgei
157  << abort(FatalError);
158  }
159  #endif
160 
161  return adjacentFace[edgei][facei];
162 }
163 
164 
165 inline Foam::edge Foam::tetCell::tetEdge(const label edgei) const
166 {
167  // Warning. Ordering of edges needs to be the same for a tetrahedron
168  // class, a tetrahedron cell shape model and a tetCell
169  static const label pt0[] = {0, 0, 0, 3, 1, 3};
170  static const label pt1[] = {1, 2, 3, 1, 2, 2};
171 
172  #ifdef FULLDEBUG
173  if (edgei < 0 || edgei >= 6)
174  {
176  << "index out of range 0 -> 5. edgei = " << edgei
177  << abort(FatalError);
178  }
179  #endif
180 
181  return Foam::edge(operator[](pt0[edgei]), operator[](pt1[edgei]));
182 }
183 
184 
185 inline Foam::edge Foam::tetCell::edge(const label edgei) const
186 {
187  return tetEdge(edgei);
188 }
189 
190 
191 inline Foam::edge Foam::tetCell::reverseEdge(const label edgei) const
192 {
193  // Reverse edge. Using a copy is cheaper than inplace flip
194  return tetEdge(edgei).reverseEdge();
195 }
196 
197 
199 {
200  return tetPointRef
201  (
202  points[operator[](0)],
203  points[operator[](1)],
204  points[operator[](2)],
205  points[operator[](3)]
206  );
207 }
208 
209 
210 // ************************************************************************* //
Foam::edge::reverseEdge
edge reverseEdge() const
Return reverse edge as copy.
Definition: edgeI.H:225
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()
Default construct, with invalid point labels (-1)
Definition: tetCellI.H:31
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetPointRef.H:47
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:128
Foam::tetCell::face
Foam::triFace face(const label facei) const
Return i-th face.
Definition: tetCellI.H:82
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::tetCell::tetEdge
Foam::edge tetEdge(const label edgei) const
Return i-th edge from tet.
Definition: tetCellI.H:165
Foam::tetCell::edgeFace
label edgeFace(const label edgei) const
Return first face adjacent to the given edge.
Definition: tetCellI.H:108
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
Foam::tetCell::edge
Foam::edge edge(const label edgei) const
Return i-th edge. Identical to tetEdge but with generic name.
Definition: tetCellI.H:185
Foam::tetCell::tet
tetPointRef tet(const UList< point > &points) const
Return the tetrahedron.
Definition: tetCellI.H:198
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::tetCell::reverseEdge
Foam::edge reverseEdge(const label edgei) const
Return i-th edge reversed.
Definition: tetCellI.H:191
triFace
face triFace(3)
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:65