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