tetIndices.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-2019 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 "tetPointRef.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  //- Point on the face, *relative to the base point*, which
93  // characterises this tet on the face.
94  label tetPti_;
95 
96 
97 public:
98 
99  // Constructors
100 
101  //- Construct null
102  tetIndices();
103 
104  //- Construct from components
105  tetIndices(label celli, label facei, label tetPtI);
106 
107 
108  //- Destructor
109  ~tetIndices();
110 
111 
112  // Member Functions
113 
114  // Access
115 
116  //- Return the cell
117  inline label cell() const;
118 
119  //- Return non-const access to the cell
120  inline label& cell();
121 
122  //- Return the face
123  inline label face() const;
124 
125  //- Return non-const access to the face
126  inline label& face();
127 
128  //- Return the characterising tetPtI
129  inline label tetPt() const;
130 
131  //- Return non-const access to the characterising tetPtI
132  inline label& tetPt();
133 
134  //- Return the indices corresponding to the tri on the face for
135  // this tet. The normal of the tri points out of the cell.
136  inline triFace faceTriIs(const polyMesh& mesh) const;
137 
138  //- Return the geometry corresponding to this tet
139  inline tetPointRef tet(const polyMesh& mesh) const;
140 
141  //- Return the geometry corresponding to the tri on the face for
142  // this tet. The normal of the tri points out of the cell.
143  inline triPointRef faceTri(const polyMesh& mesh) const;
144 
145  //- Return the geometry corresponding to the tri on the face for
146  // this tet using the old positions.
147  inline triPointRef oldFaceTri(const polyMesh& mesh) const;
148 
149 
150  // Member Operators
151 
152  inline bool operator==(const tetIndices&) const;
153  inline bool operator!=(const tetIndices&) const;
154 
155 
156  // IOstream Operators
157 
158  friend Istream& operator>>(Istream&, tetIndices&);
159  friend Ostream& operator<<(Ostream&, const tetIndices&);
160 };
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #include "tetIndicesI.H"
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 #endif
174 
175 // ************************************************************************* //
A tetrahedron primitive.
Definition: tetrahedron.H:61
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:52
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.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
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:113
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
~tetIndices()
Destructor.
Definition: tetIndices.C:53
friend Istream & operator>>(Istream &, tetIndices &)
triFace faceTriIs(const polyMesh &mesh) const
Return the indices corresponding to the tri on the face for.
Definition: tetIndicesI.H:64
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:98
dynamicFvMesh & mesh
bool operator==(const tetIndices &) const
Definition: tetIndicesI.H:146
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
label cell() const
Return the cell.
Definition: tetIndicesI.H:28
Istream & operator>>(Istream &, directionInfo &)
friend Ostream & operator<<(Ostream &, const tetIndices &)
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:40
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
bool operator!=(const tetIndices &) const
Definition: tetIndicesI.H:155
Ostream & operator<<(Ostream &, const ensightPart &)
tetIndices()
Construct null.
Definition: tetIndices.C:30
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
triPointRef oldFaceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:128
Namespace for OpenFOAM.