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 Description
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "IOstreams.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 :
37  FixedList<label, 4>(-1)
38 {}
39 
40 
42 (
43  const label a,
44  const label b,
45  const label c,
46  const label d
47 )
48 {
49  operator[](0) = a;
50  operator[](1) = b;
51  operator[](2) = c;
52  operator[](3) = d;
53 }
54 
55 
56 inline Foam::tetCell::tetCell(std::initializer_list<label> list)
57 :
58  FixedList<label, 4>(list)
59 {}
60 
61 
63 :
64  FixedList<label, 4>(list)
65 {}
66 
67 
69 (
70  const labelUList& list,
71  const FixedList<label, 4>& indices
72 )
73 :
74  FixedList<label, 4>(list, indices)
75 {}
76 
77 
79 :
80  FixedList<label, 4>(is)
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
86 inline Foam::triFace Foam::tetCell::face(const label facei) const
87 {
88  // Warning. Ordering of faces needs to be the same for a tetrahedron
89  // class, a tetrahedron cell shape model and a tetCell
90  static const label a[4] = {1, 0, 0, 0};
91  static const label b[4] = {2, 3, 1, 2};
92  static const label c[4] = {3, 2, 3, 1};
93 
94  #ifdef FULLDEBUG
95  if (facei < 0 || facei >= 4)
96  {
98  << "index out of range 0 -> 3. facei = " << facei
99  << abort(FatalError);
100  }
101  #endif
102 
103  return triFace
104  (
105  operator[](a[facei]),
106  operator[](b[facei]),
107  operator[](c[facei])
108  );
109 }
110 
111 
112 inline Foam::label Foam::tetCell::edgeFace(const label edgei) const
113 {
114  // Warning. Ordering of faces needs to be the same for a tetrahedron
115  // class, a tetrahedron cell shape model and a tetCell
116  static const label edgeFaces[6] = {2, 3, 1, 0, 0, 1};
117 
118  #ifdef FULLDEBUG
119  if (edgei < 0 || edgei >= 6)
120  {
122  << "edge index out of range 0 -> 5. edgei = " << edgei
123  << abort(FatalError);
124  }
125  #endif
126 
127  return edgeFaces[edgei];
128 }
129 
130 
131 inline Foam::label Foam::tetCell::edgeAdjacentFace
132 (
133  const label edgei,
134  const label facei
135 ) const
136 {
137  // Warning. Ordering of faces needs to be the same for a tetrahedron
138  // class, a tetrahedron cell shape model and a tetCell
139  static const label adjacentFace[6][4] =
140  {
141  {-1, -1, 3, 2},
142  {-1, 3, -1, 1},
143  {-1, 2, 1, -1},
144  {2, -1, 0, -1},
145  {3, -1, -1, 0},
146  {1, 0, -1, -1}
147  };
148 
149  #ifdef FULLDEBUG
150  if (facei < 0 || facei >= 4)
151  {
153  << "face index out of range 0 -> 3. facei = " << facei
154  << abort(FatalError);
155  }
156 
157  if (edgei < 0 || edgei >= 6)
158  {
160  << "edge index out of range 0 -> 5. edgei = " << edgei
161  << abort(FatalError);
162  }
163  #endif
164 
165  return adjacentFace[edgei][facei];
166 }
167 
168 
169 inline Foam::edge Foam::tetCell::tetEdge(const label edgei) const
170 {
171  // Warning. Ordering of edges needs to be the same for a tetrahedron
172  // class, a tetrahedron cell shape model and a tetCell
173  static const label pt0[] = {0, 0, 0, 3, 1, 3};
174  static const label pt1[] = {1, 2, 3, 1, 2, 2};
175 
176  #ifdef FULLDEBUG
177  if (edgei < 0 || edgei >= 6)
178  {
180  << "index out of range 0 -> 5. edgei = " << edgei
181  << abort(FatalError);
182  }
183  #endif
184 
185  return edge(operator[](pt0[edgei]), operator[](pt1[edgei]));
186 }
187 
188 
190 {
191  return tetPointRef
192  (
193  points[operator[](0)],
194  points[operator[](1)],
195  points[operator[](2)],
196  points[operator[](3)]
197  );
198 }
199 
200 
201 // ************************************************************************* //
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()
Default construct, with invalid point labels (-1)
Definition: tetCellI.H:35
Foam::tetPointRef
tetrahedron< point, const point & > tetPointRef
Definition: tetPointRef.H:46
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:132
Foam::tetCell::face
triFace face(const label facei) const
Return i-th face.
Definition: tetCellI.H:86
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:169
Foam::tetCell::edgeFace
label edgeFace(const label edgei) const
Return first face adjacent to the given edge.
Definition: tetCellI.H:112
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
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:189
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.
triFace
face triFace(3)
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:65