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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::tetIndices
28 
29 Description
30  Storage and named access for the indices of a tet which is part of
31  the decomposition of a cell.
32 
33  Tets are designated by
34  - cell (of course)
35  - face on cell
36  - three points on face (faceBasePt, facePtA, facePtB)
37  When constructing from a mesh and index in the face (tetPtI):
38  - faceBasePt is the mesh.tetBasePtIs() base point
39  - facePtA is tetPtI away from faceBasePt
40  - facePtB is next one after/before facePtA
41  e.g.:
42 
43  +---+
44  |2 /|
45  | / |
46  |/ 1| <- tetPt (so 1 for first triangle, 2 for second)
47  +---+
48  ^
49  faceBasePt
50 
51 SourceFiles
52  tetIndicesI.H
53  tetIndices.C
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #ifndef tetIndices_H
58 #define tetIndices_H
59 
60 #include "label.H"
61 #include "tetPointRef.H"
62 #include "triPointRef.H"
63 #include "polyMesh.H"
64 #include "triFace.H"
65 #include "face.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 // Forward declaration of friend functions and operators
73 
74 class tetIndices;
75 
76 Istream& operator>>(Istream&, tetIndices&);
77 Ostream& operator<<(Ostream&, const tetIndices&);
78 
79 
80 /*---------------------------------------------------------------------------*\
81  Class tetIndices Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 class tetIndices
85 {
86  // Private data
87 
88  //- Cell that this is a decomposed tet of
89  label celli_;
90 
91  //- Face that holds this decomposed tet
92  label facei_;
93 
94  //- Point on the face, *relative to the base point*, which
95  // characterises this tet on the face.
96  label tetPti_;
97 
98 
99  // Private static data
100 
101  //- Maximum number of bad tet warnings
102  static const label maxNWarnings;
103 
104  //- Current number of bad tet warnings. Warnings stop when this reaches
105  // the maximum number.
106  static label nWarnings;
107 
108 
109 public:
110 
111  // Constructors
112 
113  //- Construct null
114  tetIndices();
115 
116  //- Construct from components
117  tetIndices(label celli, label facei, label tetPtI);
118 
119 
120  //- Destructor
121  ~tetIndices();
122 
123 
124  // Member Functions
125 
126  // Access
127 
128  //- Return the cell
129  inline label cell() const;
130 
131  //- Return non-const access to the cell
132  inline label& cell();
133 
134  //- Return the face
135  inline label face() const;
136 
137  //- Return non-const access to the face
138  inline label& face();
139 
140  //- Return the characterising tetPtI
141  inline label tetPt() const;
142 
143  //- Return non-const access to the characterising tetPtI
144  inline label& tetPt();
145 
146  //- Return the indices corresponding to the tri on the face for
147  // this tet. The normal of the tri points out of the cell
148  inline triFace faceTriIs
149  (
150  const polyMesh& mesh,
151  const bool warn = true
152  ) const;
153 
154  //- Return the local indices corresponding to the tri on the face
155  // for this tet. The normal of the tri points out of the cell
156  inline triFace triIs
157  (
158  const polyMesh& mesh,
159  const bool warn = true
160  ) const;
161 
162  //- Return the geometry corresponding to this tet
163  inline tetPointRef tet(const polyMesh& mesh) const;
164 
165  //- Return the geometry corresponding to the tri on the face for
166  // this tet. The normal of the tri points out of the cell
167  inline triPointRef faceTri(const polyMesh& mesh) const;
168 
169  //- Return the geometry corresponding to the tri on the face for
170  // this tet using the old positions
171  inline triPointRef oldFaceTri(const polyMesh& mesh) const;
172 
173 
174  // Member Operators
175 
176  inline bool operator==(const tetIndices&) const;
177  inline bool operator!=(const tetIndices&) const;
178 
179 
180  // IOstream Operators
181 
182  friend Istream& operator>>(Istream&, tetIndices&);
183  friend Ostream& operator<<(Ostream&, const tetIndices&);
184 };
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 } // End namespace Foam
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 #include "tetIndicesI.H"
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #endif
198 
199 // ************************************************************************* //
Foam::tetIndices::tet
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:155
Foam::tetIndices::triIs
triFace triIs(const polyMesh &mesh, const bool warn=true) const
Return the local indices corresponding to the tri on the face.
Definition: tetIndicesI.H:112
Foam::tetIndices::operator==
bool operator==(const tetIndices &) const
Definition: tetIndicesI.H:203
Foam::tetIndices::cell
label cell() const
Return the cell.
Definition: tetIndicesI.H:31
Foam::tetIndices::face
label face() const
Return the face.
Definition: tetIndicesI.H:43
face.H
Foam::tetIndices::oldFaceTri
triPointRef oldFaceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:185
tetPointRef.H
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
triPointRef.H
triFace.H
Foam::tetIndices::operator<<
friend Ostream & operator<<(Ostream &, const tetIndices &)
polyMesh.H
Foam::tetIndices::tetIndices
tetIndices()
Construct null.
Definition: tetIndices.C:39
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::tetIndices::faceTri
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:170
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
tetIndicesI.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::tetIndices::faceTriIs
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
Return the indices corresponding to the tri on the face for.
Definition: tetIndicesI.H:68
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::tetIndices::tetPt
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:55
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
label.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::tetIndices::operator>>
friend Istream & operator>>(Istream &, tetIndices &)
Foam::tetIndices::operator!=
bool operator!=(const tetIndices &) const
Definition: tetIndicesI.H:212
Foam::tetIndices::~tetIndices
~tetIndices()
Destructor.
Definition: tetIndices.C:62
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:63