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-2015 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 /*---------------------------------------------------------------------------*\
71  Class tetIndices Declaration
72 \*---------------------------------------------------------------------------*/
73 
74 class tetIndices
75 {
76  // Private data
77 
78  //- Cell that this is a decomposed tet of
79  label cellI_;
80 
81  //- Face that holds this decomposed tet
82  label faceI_;
83 
84  //- Base point on the face
85  label faceBasePtI_;
86 
87  //- Point on the face such that the right-hand circulation
88  // {faceBasePtI_, facePtAI_, facePtBI_}
89  // forms a triangle that points out of the tet
90  label facePtAI_;
91 
92  //- Point on the face such that the right-hand circulation
93  // {faceBasePtI_, facePtAI_, facePtBI_}
94  // forms a triangle that points out of the tet
95  label facePtBI_;
96 
97  //- Point on the face, *relative to the base point*, which
98  // characterises this tet on the face.
99  label tetPtI_;
100 
101 
102 public:
103 
104  // Constructors
105 
106  //- Construct null
107  tetIndices();
108 
109  //- Construct from components
110  tetIndices
111  (
112  label cellI,
113  label faceI,
114  label faceBasePtI,
115  label facePtAI,
116  label facePtBI,
117  label tetPtI
118  );
119 
120  //- Construct from cell, face, pt description of tet
121  tetIndices
122  (
123  label cellI,
124  label faceI,
125  label tetPtI,
126  const polyMesh& mesh
127  );
128 
129 
130  //- Destructor
131  ~tetIndices();
132 
133 
134  // Member Functions
135 
136  // Access
137 
138  //- Return the cell
139  inline label cell() const;
140 
141  //- Return the face
142  inline label face() const;
143 
144  //- Return the face base point
145  inline label faceBasePt() const;
146 
147  //- Return face point A
148  inline label facePtA() const;
149 
150  //- Return face point B
151  inline label facePtB() const;
152 
153  //- Return the characterising tetPtI
154  inline label tetPt() const;
155 
156  //- Return the geometry corresponding to this tet from the
157  // supplied mesh
158  inline tetPointRef tet(const polyMesh& mesh) const;
159 
160  //- Return the geometry corresponding to this tet from the
161  // supplied mesh using the old positions
162  inline tetPointRef oldTet(const polyMesh& mesh) const;
163 
164  //- Return the geometry corresponding to the tri on the
165  // mesh face for this tet from the supplied mesh. Normal of
166  // the tri points out of the cell.
167  inline triPointRef faceTri(const polyMesh& mesh) const;
168 
169  //- Return the point indices corresponding to the tri on the mesh
170  // face for this tet from the supplied mesh. Normal of
171  // the tri points out of the cell.
172  inline triFace faceTriIs(const polyMesh& mesh) const;
173 
174  //- Return the geometry corresponding to the tri on the
175  // mesh face for this tet from the supplied mesh using
176  // the old position
177  inline triPointRef oldFaceTri(const polyMesh& mesh) const;
178 
179 
180  // Edit
181 
182  //- Return non-const access to the cell
183  inline label& cell();
184 
185  //- Return non-const access to the face
186  inline label& face();
187 
188  //- Return non-const access to the face base point
189  inline label& faceBasePt();
190 
191  //- Return non-const access to face point A
192  inline label& facePtA();
193 
194  //- Return non-const access to face point B
195  inline label& facePtB();
196 
197  //- Return non-const access to the characterising tetPtI
198  inline label& tetPt();
199 
200 
201  // Member Operators
202 
203  inline bool operator==(const tetIndices&) const;
204  inline bool operator!=(const tetIndices&) const;
205 
206 
207  // IOstream Operators
208 
209  friend Istream& operator>>(Istream&, tetIndices&);
210  friend Ostream& operator<<(Ostream&, const tetIndices&);
211 };
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace Foam
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #include "tetIndicesI.H"
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 #endif
225 
226 // ************************************************************************* //
label cell() const
Return the cell.
Definition: tetIndicesI.H:30
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
label faceBasePt() const
Return the face base point.
Definition: tetIndicesI.H:42
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
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
dynamicFvMesh & mesh
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:57
Namespace for OpenFOAM.
bool operator==(const tetIndices &) const
Definition: tetIndicesI.H:194
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet from the.
Definition: tetIndicesI.H:66
friend Istream & operator>>(Istream &, tetIndices &)
tetPointRef oldTet(const polyMesh &mesh) const
Return the geometry corresponding to this tet from the.
Definition: tetIndicesI.H:84
A tetrahedron primitive.
Definition: tetrahedron.H:62
bool operator!=(const tetIndices &) const
Definition: tetIndicesI.H:208
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
tetIndices()
Construct null.
Definition: tetIndices.C:30
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the.
Definition: tetIndicesI.H:109
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
label facePtA() const
Return face point A.
Definition: tetIndicesI.H:48
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:60
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:73
~tetIndices()
Destructor.
Definition: tetIndices.C:102
label facePtB() const
Return face point B.
Definition: tetIndicesI.H:54
label face() const
Return the face.
Definition: tetIndicesI.H:36