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-2021 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 
74 template<unsigned AnyNum>
76 (
77  const FixedList<label, AnyNum>& list,
78  const FixedList<label, 4>& indices
79 )
80 :
81  FixedList<label, 4>(list, indices)
82 {}
83 
84 
86 :
87  FixedList<label, 4>(is)
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
93 inline Foam::triFace Foam::tetCell::face(const label facei) const
94 {
95  #ifdef FULLDEBUG
96  if (facei < 0 || facei >= tetCell::nFaces())
97  {
99  << "Face index (" << facei << ") out of range 0..3\n"
100  << abort(FatalError);
101  }
102  #endif
103 
104  return Foam::triFace
105  (
106  (*this)[modelFaces_[facei][0]],
107  (*this)[modelFaces_[facei][1]],
108  (*this)[modelFaces_[facei][2]]
109  );
110 }
111 
112 
113 inline Foam::label Foam::tetCell::edgeFace(const label edgei) const
114 {
115  // Warning. Ordering of faces needs to be the same for a tetrahedron
116  // class, a tetrahedron cell shape model and a tetCell
117  static const label edgeFaces[6] = {2, 3, 1, 0, 0, 1};
118 
119  #ifdef FULLDEBUG
120  if (edgei < 0 || edgei >= tetCell::nEdges())
121  {
123  << "Edge index (" << edgei << ") out of range 0..5\n"
124  << abort(FatalError);
125  }
126  #endif
127 
128  return edgeFaces[edgei];
129 }
130 
131 
132 inline Foam::label Foam::tetCell::edgeAdjacentFace
133 (
134  const label edgei,
135  const label facei
136 ) const
137 {
138  // Warning. Ordering of faces needs to be the same for a tetrahedron
139  // class, a tetrahedron cell shape model and a tetCell
140  static const label adjacentFace[6][4] =
141  {
142  {-1, -1, 3, 2},
143  {-1, 3, -1, 1},
144  {-1, 2, 1, -1},
145  {2, -1, 0, -1},
146  {3, -1, -1, 0},
147  {1, 0, -1, -1}
148  };
149 
150  #ifdef FULLDEBUG
151  if (facei < 0 || facei >= tetCell::nFaces())
152  {
154  << "Face index (" << facei << ") out of range 0..3\n"
155  << abort(FatalError);
156  }
157 
158  if (edgei < 0 || edgei >= tetCell::nEdges())
159  {
161  << "Edge index (" << edgei << ") out of range 0..5\n"
162  << abort(FatalError);
163  }
164  #endif
165 
166  return adjacentFace[edgei][facei];
167 }
168 
169 
170 inline Foam::edge Foam::tetCell::edge(const label edgei) const
171 {
172  #ifdef FULLDEBUG
173  if (edgei < 0 || edgei >= tetCell::nEdges())
174  {
176  << "Edge index (" << edgei << ") out of range 0..5\n"
177  << abort(FatalError);
178  }
179  #endif
180 
181  return Foam::edge
182  (
183  (*this)[modelEdges_[edgei][0]],
184  (*this)[modelEdges_[edgei][1]]
185  );
186 }
187 
188 
189 inline Foam::edge Foam::tetCell::reverseEdge(const label edgei) const
190 {
191  // Reverse edge. Using a copy is cheaper than inplace flip
192  return this->edge(edgei).reverseEdge();
193 }
194 
195 
197 (
198  const UList<point>& meshPoints
199 ) const
200 {
201  return pointField(List<point>(meshPoints, *this));
202 }
203 
204 
206 (
207  const UList<point>& meshPoints
208 ) const
209 {
210  return tetPointRef(meshPoints, *this);
211 }
212 
213 
214 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::tetCell::nEdges
static constexpr label nEdges() noexcept
Number of edges for TET.
Definition: tetCell.H:127
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::points
pointField points(const UList< point > &meshPoints) const
The points corresponding to this shape.
Definition: tetCellI.H:197
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::tetCell::tet
tetPointRef tet(const UList< point > &meshPoints) const
Return the tetrahedron.
Definition: tetCellI.H:206
Foam::tetCell::nFaces
static constexpr label nFaces() noexcept
Number of faces for TET.
Definition: tetCell.H:133
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
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:133
Foam::tetCell::face
Foam::triFace face(const label facei) const
Return i-th face.
Definition: tetCellI.H:93
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::tetCell::edgeFace
label edgeFace(const label edgei) const
Return first face adjacent to the given edge.
Definition: tetCellI.H:113
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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.
Definition: tetCellI.H:170
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::UList< label >
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:189
triFace
face triFace(3)
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:63