cellShapeI.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-2021 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
27\*---------------------------------------------------------------------------*/
28
29#include "Istream.H"
30#include "cell.H"
31#include "cellModel.H"
32#include "IndirectList.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36inline constexpr Foam::cellShape::cellShape() noexcept
37:
38 labelList(),
39 m(nullptr)
40{}
41
42
44(
45 const cellModel& model,
46 const labelUList& labels,
47 const bool doCollapse
48)
49:
50 labelList(labels),
51 m(&model)
52{
53 if (doCollapse)
54 {
55 collapse();
56 }
57}
58
59
60template<unsigned N>
62(
63 const cellModel& model,
64 const FixedList<label, N>& labels,
65 const bool doCollapse
66)
67:
68 labelList(labels),
69 m(&model)
70{
71 if (doCollapse)
72 {
73 collapse();
74 }
75}
76
77
79(
80 const cellModel& model,
81 labelList&& labels,
82 const bool doCollapse
83)
84:
85 labelList(std::move(labels)),
86 m(&model)
87{
88 if (doCollapse)
89 {
90 collapse();
91 }
92}
93
94
96(
97 const word& modelName,
98 const labelUList& labels,
99 const bool doCollapse
100)
101:
102 labelList(labels),
103 m(cellModel::ptr(modelName))
104{
105 if (doCollapse)
106 {
107 collapse();
108 }
109}
110
111
113{
114 is >> *this;
115}
116
117
119{
120 return autoPtr<cellShape>::New(*this);
121}
122
123
124// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125
127{
128 return *m;
129}
130
131
132inline Foam::label Foam::cellShape::nPoints() const noexcept
133{
134 return labelList::size();
135}
136
137
138inline Foam::label Foam::cellShape::nEdges() const
139{
140 return m->nEdges();
141}
142
143
144inline Foam::label Foam::cellShape::nFaces() const
145{
146 return m->nFaces();
147}
148
149
151(
152 const UList<point>& meshPoints
153) const
154{
155 return pointField(UIndirectList<point>(meshPoints, *this));
156}
157
158
160(
161 const faceList& allFaces,
162 const cell& cFaces
163) const
164{
165 // Faces in model order
166 faceList localFaces(faces());
167
168 // Do linear match (usually cell shape is low complexity)
169
170 labelList modelToMesh(localFaces.size(), -1);
171
172 forAll(localFaces, i)
173 {
174 const Foam::face& localF = localFaces[i];
175
176 for (const label meshFacei : cFaces)
177 {
178 if (allFaces[meshFacei] == localF)
179 {
180 modelToMesh[i] = meshFacei;
181 break;
182 }
183 }
184 }
185
186 return modelToMesh;
187}
188
189
191(
192 const edgeList& allEdges,
193 const labelList& cEdges
194) const
195{
196 // Edges in model order
197 edgeList localEdges(edges());
198
199 // Do linear match (usually cell shape is low complexity)
200
201 labelList modelToMesh(localEdges.size(), -1);
202
203 forAll(localEdges, i)
204 {
205 const Foam::edge& e = localEdges[i];
206
207 for (const label edgei : cEdges)
208 {
209 if (allEdges[edgei] == e)
210 {
211 modelToMesh[i] = edgei;
212 break;
213 }
214 }
215 }
216
217 return modelToMesh;
218}
219
220
221inline Foam::face Foam::cellShape::face(const label modelFacei) const
222{
223 return m->face(modelFacei, *this);
224}
225
226
228{
229 return m->faces(*this);
230}
231
232
234{
235 const faceList oldFaces(faces());
236
237 faceList newFaces(oldFaces.size());
238 label newFacei = 0;
239
240 for (const Foam::face& f : oldFaces)
241 {
242 Foam::face& newF = newFaces[newFacei];
243 newF.resize(f.size());
244
245 label newFp = 0;
246 label prevVerti = -1;
247
248 for (const label verti : f)
249 {
250 if (verti != prevVerti)
251 {
252 newF[newFp++] = verti;
253 prevVerti = verti;
254 }
255 }
256
257 if ((newFp > 1) && (newF[newFp-1] == newF[0]))
258 {
259 --newFp;
260 }
261
262 if (newFp > 2)
263 {
264 // Size face and go to next one
265 newF.resize(newFp);
266 ++newFacei;
267 }
268 }
269 newFaces.resize(newFacei);
270
271 return newFaces;
272}
273
274
275inline Foam::edge Foam::cellShape::edge(const label modelEdgei) const
276{
277 return m->edge(modelEdgei, *this);
278}
279
280
282{
283 return m->edges(*this);
284}
285
286
288{
289 return m->centre(*this, points);
290}
291
292
293inline Foam::scalar Foam::cellShape::mag(const UList<point>& points) const
294{
295 return m->mag(*this, points);
296}
297
298
300(
301 const cellModel& model,
302 const labelUList& labels,
303 const bool doCollapse
304)
305{
306 static_cast<labelList&>(*this) = labels;
307 m = &model;
308
309 if (doCollapse)
310 {
311 collapse();
312 }
313}
314
315
316template<unsigned N>
318(
319 const cellModel& model,
320 const FixedList<label, N>& labels,
321 const bool doCollapse
322)
323{
324 static_cast<labelList&>(*this) = labels;
325 m = &model;
326
327 if (doCollapse)
328 {
329 collapse();
330 }
331}
332
333
334// ************************************************************************* //
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
scalar centre() const
Mid-point location, zero for an empty list.
Definition: PDRblockI.H:67
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:73
const cellModel & model() const
Model reference.
Definition: cellShapeI.H:126
autoPtr< cellShape > clone() const
Clone.
Definition: cellShapeI.H:118
void collapse()
Collapse shape to correct one after removing duplicate vertices.
Definition: cellShape.C:33
labelList meshFaces(const faceList &allFaces, const cell &cFaces) const
Mesh face labels of this cell (in order of model)
Definition: cellShapeI.H:160
label nPoints() const noexcept
Number of points.
Definition: cellShapeI.H:132
faceList faces() const
Faces of this cell.
Definition: cellShapeI.H:227
constexpr cellShape() noexcept
Default construct. Empty shape, no cell model.
Definition: cellShapeI.H:36
edgeList edges() const
Edges of this shape.
Definition: cellShapeI.H:281
faceList collapsedFaces() const
Collapsed faces of this cell.
Definition: cellShapeI.H:233
label nEdges() const
Number of edges.
Definition: cellShapeI.H:138
label nFaces() const
Number of faces.
Definition: cellShapeI.H:144
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
edge()
Default construct, with invalid point labels (-1)
Definition: edgeI.H:41
void reset()
Reset to defaults.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
constexpr face() noexcept=default
Default construct.
Computes the magnitude of an input field.
Definition: mag.H:153
mag(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: mag.C:62
const labelList & meshEdges() const
Return global edge index for local edges.
Definition: polyPatch.C:385
A class for handling words, derived from Foam::string.
Definition: word.H:68
const pointField & points
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
const direction noexcept
Definition: Scalar.H:223
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333