Loading...
Searching...
No Matches
cellShape.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-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::cellShape
29
30Description
31 An analytical geometric cellShape.
32
33 The optional collapse functionality changes the cellModel to the
34 correct type after removing any duplicate points.
35
36SourceFiles
37 cellShapeI.H
38 cellShape.C
39 cellShapeIO.C
40 cellShapeEqual.C
41
42\*---------------------------------------------------------------------------*/
43
44#ifndef Foam_cellShape_H
45#define Foam_cellShape_H
46
47#include "pointField.H"
48#include "labelList.H"
49#include "cellModel.H"
50#include "autoPtr.H"
51#include "InfoProxy.H"
52
53// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54
55namespace Foam
56{
57
58// Forward Declarations
59class cell;
60class cellShape;
61bool operator==(const cellShape& a, const cellShape& b);
64
65template<>
67
69/*---------------------------------------------------------------------------*\
70 Class cellShape Declaration
71\*---------------------------------------------------------------------------*/
72
73class cellShape
74:
75 public labelList
76{
77 // Private Data
78
79 //- Access to the cellShape's model
80 const cellModel *m;
81
82
83public:
84
85 // Constructors
86
87 //- Default construct. Empty shape, no cell model.
88 inline constexpr cellShape() noexcept;
89
90 //- Copy construct from components
91 inline cellShape
92 (
93 const cellModel& model,
94 const labelUList& labels,
95 const bool doCollapse = false
96 );
97
98 //- Copy construct from components
99 template<unsigned N>
100 inline cellShape
101 (
102 const cellModel& model,
103 const FixedList<label, N>& labels,
104 const bool doCollapse = false
105 );
106
107 //- Move construct from components
108 inline cellShape
109 (
110 const cellModel& model,
111 labelList&& labels,
112 const bool doCollapse = false
113 );
114
115 //- Copy construct from components, lookup cellModel by type
116 inline cellShape
117 (
119 const labelUList& labels,
120 const bool doCollapse = false
121 );
122
123 //- Move construct from components, lookup cellModel by type
124 inline cellShape
125 (
127 labelList&& labels,
128 const bool doCollapse = false
129 );
130
131 //- Copy construct from components, lookup cellModel by name
132 inline cellShape
133 (
134 const word& modelName,
135 const labelUList& labels,
136 const bool doCollapse = false
137 );
138
139 //- Construct from Istream
140 inline explicit cellShape(Istream& is);
141
142 //- Clone
143 inline autoPtr<cellShape> clone() const;
144
145
146 // Member Functions
147
148 //- Model reference
149 inline const cellModel& model() const;
150
151 //- Number of points
152 inline label nPoints() const noexcept;
153
154 //- Number of edges
155 inline label nEdges() const;
156
157 //- Number of faces
158 inline label nFaces() const;
159
160 //- The points corresponding to this shape
161 inline pointField points(const UList<point>& meshPoints) const;
162
163 //- Mesh face labels of this cell (in order of model)
164 inline labelList meshFaces
165 (
166 const faceList& allFaces,
167 const cell& cFaces
168 ) const;
169
170 //- Mesh edge labels of this cell (in order of model)
171 inline labelList meshEdges
172 (
173 const edgeList& allEdges,
174 const labelList& cEdges
175 ) const;
176
177 //- The face for the specified model face
178 inline Foam::face face(const label modelFacei) const;
179
180 //- Faces of this cell
181 inline faceList faces() const;
182
183 //- Collapsed faces of this cell
184 inline faceList collapsedFaces() const;
185
186 //- The edge for the specified model edge
187 inline Foam::edge edge(const label modelEdgei) const;
188
189 //- Edges of this shape
190 inline edgeList edges() const;
191
192 //- Centroid of the cell
193 inline point centre(const UList<point>& points) const;
194
195 //- Scalar magnitude
196 inline scalar mag(const UList<point>& points) const;
197
198 //- Return the smallest vertex label used by the shape
199 //- or -1 if the shape is empty.
200 inline label min() const;
201
202 //- Return the largest vertex label used by the shape
203 //- or -1 if the shape is empty.
204 inline label max() const;
205
206 //- Reset from components
207 inline void reset
208 (
209 const cellModel& model,
210 const labelUList& labels,
211 const bool doCollapse = false
212 );
213
214 //- Reset from components
215 template<unsigned N>
216 inline void reset
217 (
218 const cellModel& model,
219 const FixedList<label, N>& labels,
220 const bool doCollapse = false
221 );
222
223 //- Collapse shape to correct one after removing duplicate vertices
224 void collapse();
225
226 //- Return info proxy,
227 //- used to print information to a stream
228 InfoProxy<cellShape> info() const noexcept { return *this; }
229
230
231 // Comparison
232
233 //- Compare cellShape vertices
234 // \return
235 // - 0: different
236 // - +1: identical values and order used
237 // - -1: identical values, but in different order
238 static int compare(const cellShape& a, const cellShape& b);
239
240
241 // Member Operators
242
243 //- Increment (offset) vertices by given amount
244 inline void operator+=(const label vertexOffset);
245
246
247 // Friend Operators
248
249 friend bool operator==(const cellShape& a, const cellShape& b);
250
251
252 // IOstream Operators
253
254 friend Istream& operator>>(Istream& is, cellShape& s);
255 friend Ostream& operator<<(Ostream& os, const cellShape& s);
256};
257
258
259// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260
261} // End namespace Foam
262
263// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264
265#include "cellShapeI.H"
266
267// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268
269#endif
270
271// ************************************************************************* //
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
A helper class for outputting values to Ostream.
Definition InfoProxy.H:49
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Maps a geometry to a set of cell primitives.
Definition cellModel.H:73
modelType
Enumeration of commonly used cellModel types.
Definition cellModel.H:82
An analytical geometric cellShape.
Definition cellShape.H:71
const cellModel & model() const
Model reference.
Definition cellShapeI.H:153
labelList meshEdges(const edgeList &allEdges, const labelList &cEdges) const
Mesh edge labels of this cell (in order of model).
Definition cellShapeI.H:218
label min() const
Return the smallest vertex label used by the shape or -1 if the shape is empty.
Definition cellShapeI.H:326
friend bool operator==(const cellShape &a, const cellShape &b)
scalar mag(const UList< point > &points) const
Scalar magnitude.
Definition cellShapeI.H:320
void reset(const cellModel &model, const labelUList &labels, const bool doCollapse=false)
Reset from components.
Definition cellShapeI.H:359
autoPtr< cellShape > clone() const
Clone.
Definition cellShapeI.H:145
void collapse()
Collapse shape to correct one after removing duplicate vertices.
Definition cellShape.C:26
Foam::face face(const label modelFacei) const
The face for the specified model face.
Definition cellShapeI.H:248
pointField points(const UList< point > &meshPoints) const
The points corresponding to this shape.
Definition cellShapeI.H:178
labelList meshFaces(const faceList &allFaces, const cell &cFaces) const
Mesh face labels of this cell (in order of model).
Definition cellShapeI.H:187
void operator+=(const label vertexOffset)
Increment (offset) vertices by given amount.
Definition cellShapeI.H:395
point centre(const UList< point > &points) const
Centroid of the cell.
Definition cellShapeI.H:314
friend Ostream & operator<<(Ostream &os, const cellShape &s)
static int compare(const cellShape &a, const cellShape &b)
Compare cellShape vertices.
label nPoints() const noexcept
Number of points.
Definition cellShapeI.H:159
faceList faces() const
Faces of this cell.
Definition cellShapeI.H:254
constexpr cellShape() noexcept
Default construct. Empty shape, no cell model.
Definition cellShapeI.H:29
edgeList edges() const
Edges of this shape.
Definition cellShapeI.H:308
Foam::edge edge(const label modelEdgei) const
The edge for the specified model edge.
Definition cellShapeI.H:302
label max() const
Return the largest vertex label used by the shape or -1 if the shape is empty.
Definition cellShapeI.H:342
faceList collapsedFaces() const
Collapsed faces of this cell.
Definition cellShapeI.H:260
label nEdges() const
Number of edges.
Definition cellShapeI.H:165
friend Istream & operator>>(Istream &is, cellShape &s)
label nFaces() const
Number of faces.
Definition cellShapeI.H:171
InfoProxy< cellShape > info() const noexcept
Return info proxy, used to print information to a stream.
Definition cellShape.H:283
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling words, derived from Foam::string.
Definition word.H:66
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< label > labelList
A List of labels.
Definition List.H:62
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Istream & operator>>(Istream &, directionInfo &)
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
vectorField pointField
pointField is a vectorField.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
volScalarField & b