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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::cellShape
29 
30 Description
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 
36 SourceFiles
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 
55 namespace Foam
56 {
57 
58 // Forward Declarations
59 class cell;
60 class cellShape;
61 bool operator==(const cellShape& a, const cellShape& b);
62 Istream& operator>>(Istream& is, cellShape& s);
63 Ostream& operator<<(Ostream& os, const cellShape& s);
64 
65 template<>
66 Ostream& operator<<(Ostream& os, const InfoProxy<cellShape>& ip);
67 
68 
69 /*---------------------------------------------------------------------------*\
70  Class cellShape Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 class cellShape
74 :
75  public labelList
76 {
77  // Private Data
78 
79  //- Access to the cellShape's model
80  const cellModel *m;
81 
82 
83 public:
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  //- Reset from components
199  inline void reset
200  (
201  const cellModel& model,
202  const labelUList& labels,
203  const bool doCollapse = false
204  );
205 
206  //- Reset from components
207  template<unsigned N>
208  inline void reset
209  (
210  const cellModel& model,
211  const FixedList<label, N>& labels,
212  const bool doCollapse = false
213  );
214 
215  //- Collapse shape to correct one after removing duplicate vertices
216  void collapse();
217 
218  //- Return info proxy,
219  //- used to print information to a stream
220  InfoProxy<cellShape> info() const noexcept { return *this; }
221 
222 
223  // Friend Operators
224 
225  friend bool operator==(const cellShape& a, const cellShape& b);
226 
227 
228  // IOstream Operators
229 
230  friend Istream& operator>>(Istream& is, cellShape& s);
231  friend Ostream& operator<<(Ostream& os, const cellShape& s);
232 };
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace Foam
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #include "cellShapeI.H"
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 #endif
246 
247 // ************************************************************************* //
label nPoints() const noexcept
Number of points.
Definition: cellShapeI.H:159
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
void collapse()
Collapse shape to correct one after removing duplicate vertices.
Definition: cellShape.C:26
An analytical geometric cellShape.
Definition: cellShape.H:68
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
InfoProxy< cellShape > info() const noexcept
Return info proxy, used to print information to a stream.
Definition: cellShape.H:271
autoPtr< cellShape > clone() const
Clone.
Definition: cellShapeI.H:145
labelList meshEdges(const edgeList &allEdges, const labelList &cEdges) const
Mesh edge labels of this cell (in order of model)
Definition: cellShapeI.H:218
friend Istream & operator>>(Istream &is, cellShape &s)
friend bool operator==(const cellShape &a, const cellShape &b)
friend Ostream & operator<<(Ostream &os, const cellShape &s)
faceList faces() const
Faces of this cell.
Definition: cellShapeI.H:254
label nEdges() const
Number of edges.
Definition: cellShapeI.H:165
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:59
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
labelList meshFaces(const faceList &allFaces, const cell &cFaces) const
Mesh face labels of this cell (in order of model)
Definition: cellShapeI.H:187
A class for handling words, derived from Foam::string.
Definition: word.H:63
Istream & operator>>(Istream &, directionInfo &)
void reset(const cellModel &model, const labelUList &labels, const bool doCollapse=false)
Reset from components.
Definition: cellShapeI.H:327
pointField points(const UList< point > &meshPoints) const
The points corresponding to this shape.
Definition: cellShapeI.H:178
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
Foam::edge edge(const label modelEdgei) const
The edge for the specified model edge.
Definition: cellShapeI.H:302
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
faceList collapsedFaces() const
Collapsed faces of this cell.
Definition: cellShapeI.H:260
OBJstream os(runTime.globalPath()/outputName)
edgeList edges() const
Edges of this shape.
Definition: cellShapeI.H:308
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
modelType
Enumeration of commonly used cellModel types.
Definition: cellModel.H:81
const cellModel & model() const
Model reference.
Definition: cellShapeI.H:153
A helper class for outputting values to Ostream.
Definition: ensightCells.H:43
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:72
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
label nFaces() const
Number of faces.
Definition: cellShapeI.H:171
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
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))
constexpr cellShape() noexcept
Default construct. Empty shape, no cell model.
Definition: cellShapeI.H:29
scalar mag(const UList< point > &points) const
Scalar magnitude.
Definition: cellShapeI.H:320
Namespace for OpenFOAM.
Foam::face face(const label modelFacei) const
The face for the specified model face.
Definition: cellShapeI.H:248
point centre(const UList< point > &points) const
Centroid of the cell.
Definition: cellShapeI.H:314