tetIndices.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-2017 OpenFOAM Foundation
9  Copyright (C) 2022 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 Class
28  Foam::tetIndices
29 
30 Description
31  Storage and named access for the indices of a tet which is part of
32  the decomposition of a cell.
33 
34  Tets are designated by
35  - cell (of course)
36  - face on cell
37  - three points on face (faceBasePt, facePtA, facePtB)
38  When constructing from a mesh and index in the face (tetPtI):
39  - faceBasePt is the mesh.tetBasePtIs() base point
40  - facePtA is tetPtI away from faceBasePt
41  - facePtB is next one after/before facePtA
42  e.g.:
43 
44  +---+
45  |2 /|
46  | / |
47  |/ 1| <- tetPt (so 1 for first triangle, 2 for second)
48  +---+
49  ^
50  faceBasePt
51 
52 SourceFiles
53  tetIndicesI.H
54  tetIndices.C
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #ifndef Foam_tetIndices_H
59 #define Foam_tetIndices_H
60 
61 #include "label.H"
62 #include "tetrahedron.H"
63 #include "triangle.H"
64 #include "polyMesh.H"
65 #include "triFace.H"
66 #include "face.H"
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 
73 // Forward Declarations
74 class tetIndices;
75 
76 Istream& operator>>(Istream&, tetIndices&);
77 Ostream& operator<<(Ostream&, const tetIndices&);
78 
79 /*---------------------------------------------------------------------------*\
80  Class tetIndices Declaration
81 \*---------------------------------------------------------------------------*/
82 
83 class tetIndices
84 {
85  // Private Data
86 
87  //- Cell that this is a decomposed tet of
88  label celli_;
89 
90  //- Face that holds this decomposed tet
91  label facei_;
92 
93  //- Point on the face, *relative to the base point*,
94  //- which characterises this tet on the face.
95  label tetPti_;
96 
97 
98  // Private Static Data
99 
100  //- Maximum number of bad tet warnings
101  static constexpr int maxNWarnings = 100;
102 
103  //- Current number of bad tet warnings.
104  //- Warnings stop when this reaches the maximum number.
105  static int nWarnings_;
106 
107 
108 public:
109 
110  // Constructors
111 
112  //- Default construct, with invalid labels (-1)
113  inline constexpr tetIndices() noexcept;
114 
115  //- Construct from components
116  inline constexpr tetIndices
117  (
118  label celli,
119  label facei,
120  label tetPointi
121  ) noexcept;
122 
123 
124  //- Destructor
125  ~tetIndices() = default;
126 
127 
128  // Member Functions
129 
130  // Access
131 
132  //- Return the cell index
133  label cell() const noexcept { return celli_; }
134 
135  //- Non-const access to the cell index
136  label& cell() noexcept { return celli_; }
137 
138  //- Return the face index
139  label face() const noexcept { return facei_; }
140 
141  //- Non-const access to the face index
142  label& face() noexcept { return facei_; }
143 
144  //- Return the characterising tet point index
145  label tetPt() const noexcept { return tetPti_; }
147  //- Non-const access to the characterising tet point index
148  label& tetPt() noexcept { return tetPti_; }
149 
150 
151  // Mesh-related
152 
153  //- The indices corresponding to the tri on the face for this tet.
154  //- The normal of the tri points out of the cell.
155  inline triFace faceTriIs
156  (
157  const polyMesh& mesh,
158  const bool warn = true
159  ) const;
160 
161  //- Local indices corresponding to the tri on the face for this tet.
162  //- The normal of the tri points out of the cell.
163  inline triFace triIs
164  (
165  const polyMesh& mesh,
166  const bool warn = true
167  ) const;
168 
169 
170  //- The tet geometry for this tet,
171  //- where point0 is the cell centre.
172  inline tetPointRef tet(const polyMesh& mesh) const;
173 
174  //- The tet geometry for this tet (using old positions),
175  //- where point0 is the cell centre.
176  inline tetPointRef oldTet(const polyMesh& mesh) const;
177 
178  //- The triangle geometry for the face for this tet.
179  //- The normal of the tri points out of the cell
180  inline triPointRef faceTri(const polyMesh& mesh) const;
181 
182  //- The triangle geometry for the face for this tet
183  //- (using old positions)
184  inline triPointRef oldFaceTri(const polyMesh& mesh) const;
185 
186  //- The x/y/z position for given barycentric coordinates
187  //- (where point0 is the cell centre).
189  (
190  const polyMesh& mesh,
191  const barycentric& bary
192  ) const;
193 
194 
195  // Other
196 
197  //- Compare tetIndices for equality.
198  //- Compares cell, face, tetPt elements in order, stopping at the
199  //- first inequality.
200  //
201  // \returns negative/zero/positive from the last element compared
202  static inline label compare
203  (
204  const tetIndices& a,
205  const tetIndices& b
206  ) noexcept;
207 
208 
209  // IOstream Operators
210 
211  friend Istream& operator>>(Istream&, tetIndices&);
212  friend Ostream& operator<<(Ostream&, const tetIndices&);
213 };
214 
215 
216 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
217 
218 //- Contiguous data for tetIndices
219 template<> struct is_contiguous<tetIndices> : std::true_type {};
220 
221 //- Contiguous label data for tetIndices
222 template<> struct is_contiguous_label<tetIndices> : std::true_type {};
223 
224 
225 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
226 
227 inline bool operator==(const tetIndices& a, const tetIndices& b) noexcept;
228 inline bool operator!=(const tetIndices& a, const tetIndices& b) noexcept;
229 
230 
231 inline bool operator<(const tetIndices& a, const tetIndices& b) noexcept
232 {
233  return (tetIndices::compare(a, b) < 0);
234 }
235 
236 inline bool operator<=(const tetIndices& a, const tetIndices& b) noexcept
237 {
238  return (tetIndices::compare(a, b) <= 0);
239 }
240 
241 inline bool operator>(const tetIndices& a, const tetIndices& b) noexcept
242 {
243  return (tetIndices::compare(a, b) > 0);
244 }
245 
246 inline bool operator>=(const tetIndices& a, const tetIndices& b) noexcept
247 {
248  return (tetIndices::compare(a, b) >= 0);
249 }
250 
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 } // End namespace Foam
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #include "tetIndicesI.H"
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 #endif
263 
264 // ************************************************************************* //
A tetrahedron primitive.
Definition: tetrahedron.H:58
~tetIndices()=default
Destructor.
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
point barycentricToPoint(const polyMesh &mesh, const barycentric &bary) const
The x/y/z position for given barycentric coordinates (where point0 is the cell centre).
Definition: tetIndicesI.H:160
triPointRef faceTri(const polyMesh &mesh) const
The triangle geometry for the face for this tet. The normal of the tri points out of the cell...
Definition: tetIndicesI.H:182
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
bool operator>(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A newer than B.
friend Istream & operator>>(Istream &, tetIndices &)
bool operator<(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A older than B.
tetPointRef tet(const polyMesh &mesh) const
The tet geometry for this tet, where point0 is the cell centre.
Definition: tetIndicesI.H:129
bool operator>=(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A same or newer than B.
bool operator<=(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A same or older than B.
dynamicFvMesh & mesh
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
label cell() const noexcept
Return the cell index.
Definition: tetIndices.H:146
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:65
Istream & operator>>(Istream &, directionInfo &)
friend Ostream & operator<<(Ostream &, const tetIndices &)
label tetPt() const noexcept
Return the characterising tet point index.
Definition: tetIndices.H:166
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:78
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
constexpr tetIndices() noexcept
Default construct, with invalid labels (-1)
Definition: tetIndicesI.H:24
tetPointRef oldTet(const polyMesh &mesh) const
The tet geometry for this tet (using old positions), where point0 is the cell centre.
Definition: tetIndicesI.H:144
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
triFace triIs(const polyMesh &mesh, const bool warn=true) const
Local indices corresponding to the tri on the face for this tet. The normal of the tri points out of ...
Definition: tetIndicesI.H:89
label face() const noexcept
Return the face index.
Definition: tetIndices.H:156
A template class to specify that a data type can be considered as being contiguous in memory...
Definition: contiguous.H:70
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:297
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
triPointRef oldFaceTri(const polyMesh &mesh) const
The triangle geometry for the face for this tet (using old positions)
Definition: tetIndicesI.H:189
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
The indices corresponding to the tri on the face for this tet. The normal of the tri points out of th...
Definition: tetIndicesI.H:48
Namespace for OpenFOAM.
static label compare(const tetIndices &a, const tetIndices &b) noexcept
Compare tetIndices for equality. Compares cell, face, tetPt elements in order, stopping at the first ...
Definition: tetIndicesI.H:198