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 | Copyright (C) 2011-2016 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 Class
25  Foam::tetIndices
26 
27 Description
28  Storage and named access for the indices of a tet which is part of
29  the decomposition of a cell.
30 
31  Tets are designated by
32  - cell (of course)
33  - face on cell
34  - three points on face (faceBasePt, facePtA, facePtB)
35  When constructing from a mesh and index in the face (tetPtI):
36  - faceBasePt is the mesh.tetBasePtIs() base point
37  - facePtA is tetPtI away from faceBasePt
38  - facePtB is next one after/before facePtA
39  e.g.:
40 
41  +---+
42  |2 /|
43  | / |
44  |/ 1| <- tetPt (so 1 for first triangle, 2 for second)
45  +---+
46  ^
47  faceBasePt
48 
49 SourceFiles
50  tetIndicesI.H
51  tetIndices.C
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef tetIndices_H
56 #define tetIndices_H
57 
58 #include "label.H"
59 #include "tetrahedron.H"
60 #include "triPointRef.H"
61 #include "polyMesh.H"
62 #include "triFace.H"
63 #include "face.H"
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 
70 // Forward declaration of friend functions and operators
71 
72 class tetIndices;
73 
74 Istream& operator>>(Istream&, tetIndices&);
75 Ostream& operator<<(Ostream&, const tetIndices&);
76 
77 
78 /*---------------------------------------------------------------------------*\
79  Class tetIndices Declaration
80 \*---------------------------------------------------------------------------*/
81 
82 class tetIndices
83 {
84  // Private data
85 
86  //- Cell that this is a decomposed tet of
87  label celli_;
88 
89  //- Face that holds this decomposed tet
90  label facei_;
91 
92  //- Base point on the face
93  label faceBasePtI_;
94 
95  //- Point on the face such that the right-hand circulation
96  // {faceBasePtI_, facePtAI_, facePtBI_}
97  // forms a triangle that points out of the tet
98  label facePtAI_;
99 
100  //- Point on the face such that the right-hand circulation
101  // {faceBasePtI_, facePtAI_, facePtBI_}
102  // forms a triangle that points out of the tet
103  label facePtBI_;
104 
105  //- Point on the face, *relative to the base point*, which
106  // characterises this tet on the face.
107  label tetPtI_;
108 
109 
110 public:
111 
112  // Constructors
113 
114  //- Construct null
115  tetIndices();
116 
117  //- Construct from components
118  tetIndices
119  (
120  label celli,
121  label facei,
122  label faceBasePtI,
123  label facePtAI,
124  label facePtBI,
125  label tetPtI
126  );
127 
128  //- Construct from cell, face, pt description of tet
129  tetIndices
130  (
131  label celli,
132  label facei,
133  label tetPtI,
134  const polyMesh& mesh
135  );
136 
137 
138  //- Destructor
139  ~tetIndices();
140 
141 
142  // Member Functions
143 
144  // Access
145 
146  //- Return the cell
147  inline label cell() const;
148 
149  //- Return the face
150  inline label face() const;
151 
152  //- Return the face base point
153  inline label faceBasePt() const;
154 
155  //- Return face point A
156  inline label facePtA() const;
157 
158  //- Return face point B
159  inline label facePtB() const;
160 
161  //- Return the characterising tetPtI
162  inline label tetPt() const;
163 
164  //- Return the geometry corresponding to this tet from the
165  // supplied mesh
166  inline tetPointRef tet(const polyMesh& mesh) const;
167 
168  //- Return the geometry corresponding to this tet from the
169  // supplied mesh using the old positions
170  inline tetPointRef oldTet(const polyMesh& mesh) const;
171 
172  //- Return the geometry corresponding to the tri on the
173  // mesh face for this tet from the supplied mesh. Normal of
174  // the tri points out of the cell.
175  inline triPointRef faceTri(const polyMesh& mesh) const;
176 
177  //- Return the point indices corresponding to the tri on the mesh
178  // face for this tet from the supplied mesh. Normal of
179  // the tri points out of the cell.
180  inline triFace faceTriIs(const polyMesh& mesh) const;
181 
182  //- Return the geometry corresponding to the tri on the
183  // mesh face for this tet from the supplied mesh using
184  // the old position
185  inline triPointRef oldFaceTri(const polyMesh& mesh) const;
186 
187 
188  // Edit
189 
190  //- Return non-const access to the cell
191  inline label& cell();
192 
193  //- Return non-const access to the face
194  inline label& face();
195 
196  //- Return non-const access to the face base point
197  inline label& faceBasePt();
198 
199  //- Return non-const access to face point A
200  inline label& facePtA();
201 
202  //- Return non-const access to face point B
203  inline label& facePtB();
204 
205  //- Return non-const access to the characterising tetPtI
206  inline label& tetPt();
207 
208 
209  // Member Operators
210 
211  inline bool operator==(const tetIndices&) const;
212  inline bool operator!=(const tetIndices&) const;
213 
214 
215  // IOstream Operators
216 
217  friend Istream& operator>>(Istream&, tetIndices&);
218  friend Ostream& operator<<(Ostream&, const tetIndices&);
219 };
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace Foam
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #include "tetIndicesI.H"
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 #endif
233 
234 // ************************************************************************* //
A tetrahedron primitive.
Definition: tetrahedron.H:62
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:57
label cell() const
Return the cell.
Definition: tetIndicesI.H:30
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
label facePtA() const
Return face point A.
Definition: tetIndicesI.H:48
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the.
Definition: tetIndicesI.H:109
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
bool operator!=(const tetIndices &) const
Definition: tetIndicesI.H:208
~tetIndices()
Destructor.
Definition: tetIndices.C:102
friend Istream & operator>>(Istream &, tetIndices &)
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:60
bool operator==(const tetIndices &) const
Definition: tetIndicesI.H:194
label faceBasePt() const
Return the face base point.
Definition: tetIndicesI.H:42
dynamicFvMesh & mesh
label facePtB() const
Return face point B.
Definition: tetIndicesI.H:54
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
Istream & operator>>(Istream &, directionInfo &)
friend Ostream & operator<<(Ostream &, const tetIndices &)
triPointRef oldFaceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the.
Definition: tetIndicesI.H:140
triFace faceTriIs(const polyMesh &mesh) const
Return the point indices corresponding to the tri on the mesh.
Definition: tetIndicesI.H:125
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
label face() const
Return the face.
Definition: tetIndicesI.H:36
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
tetPointRef oldTet(const polyMesh &mesh) const
Return the geometry corresponding to this tet from the.
Definition: tetIndicesI.H:84
Ostream & operator<<(Ostream &, const ensightPart &)
tetIndices()
Construct null.
Definition: tetIndices.C:30
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet from the.
Definition: tetIndicesI.H:66
Namespace for OpenFOAM.