Loading...
Searching...
No Matches
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-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
27\*---------------------------------------------------------------------------*/
28
29#include "Istream.H"
30#include "cell.H"
31#include "cellModel.H"
32#include "IndirectList.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
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 cellModel::modelType model,
98 const labelUList& labels,
99 const bool doCollapse
100)
101:
102 labelList(labels),
103 m(cellModel::ptr(model))
104{
105 if (doCollapse)
106 {
107 collapse();
108 }
109}
110
111
113(
114 const cellModel::modelType model,
115 labelList&& labels,
116 const bool doCollapse
117)
118:
119 labelList(std::move(labels)),
120 m(cellModel::ptr(model))
121{
122 if (doCollapse)
123 {
124 collapse();
125 }
126}
127
128
130(
131 const word& modelName,
132 const labelUList& labels,
133 const bool doCollapse
134)
135:
136 labelList(labels),
137 m(cellModel::ptr(modelName))
138{
139 if (doCollapse)
140 {
141 collapse();
142 }
143}
144
147{
148 is >> *this;
149}
150
151
154 return autoPtr<cellShape>::New(*this);
155}
156
157
158// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160inline const Foam::cellModel& Foam::cellShape::model() const
161{
162 return *m;
163}
164
166inline Foam::label Foam::cellShape::nPoints() const noexcept
167{
168 return labelList::size();
169}
170
172inline Foam::label Foam::cellShape::nEdges() const
173{
174 return m->nEdges();
175}
176
178inline Foam::label Foam::cellShape::nFaces() const
179{
180 return m->nFaces();
181}
182
183
185(
186 const UList<point>& meshPoints
187) const
188{
189 return pointField(UIndirectList<point>(meshPoints, *this));
190}
191
192
194(
195 const faceList& allFaces,
196 const cell& cFaces
197) const
198{
199 // Faces in model order
200 faceList localFaces(faces());
201
202 // Do linear match (usually cell shape is low complexity)
203
204 labelList modelToMesh(localFaces.size(), -1);
205
206 forAll(localFaces, i)
207 {
208 const auto& localF = localFaces[i];
209
210 for (const label meshFacei : cFaces)
211 {
212 if (allFaces[meshFacei] == localF)
213 {
214 modelToMesh[i] = meshFacei;
215 break;
216 }
218 }
219
220 return modelToMesh;
221}
222
223
225(
226 const edgeList& allEdges,
227 const labelList& cEdges
228) const
229{
230 // Edges in model order
231 edgeList localEdges(edges());
232
233 // Do linear match (usually cell shape is low complexity)
234
235 labelList modelToMesh(localEdges.size(), -1);
236
237 forAll(localEdges, i)
238 {
239 const auto& e = localEdges[i];
240
241 for (const label edgei : cEdges)
242 {
243 if (allEdges[edgei] == e)
244 {
245 modelToMesh[i] = edgei;
246 break;
247 }
249 }
250
251 return modelToMesh;
252}
253
255inline Foam::face Foam::cellShape::face(const label modelFacei) const
256{
257 return m->face(modelFacei, *this);
258}
259
262{
263 return m->faces(*this);
264}
265
266
268{
269 const faceList oldFaces(faces());
270
271 faceList newFaces(oldFaces.size());
272 label newFacei = 0;
273
274 for (const auto& f : oldFaces)
275 {
276 auto& newF = newFaces[newFacei];
277 newF.resize(f.size());
278
279 label newFp = 0;
280 label prevVerti = -1;
281
282 for (const label verti : f)
283 {
284 if (verti != prevVerti)
285 {
286 newF[newFp++] = verti;
287 prevVerti = verti;
288 }
289 }
290
291 if ((newFp > 1) && (newF[newFp-1] == newF[0]))
292 {
293 --newFp;
294 }
295
296 if (newFp > 2)
297 {
298 // Size face and go to next one
299 newF.resize(newFp);
300 ++newFacei;
301 }
303 newFaces.resize(newFacei);
304
305 return newFaces;
306}
307
309inline Foam::edge Foam::cellShape::edge(const label modelEdgei) const
310{
311 return m->edge(modelEdgei, *this);
312}
313
316{
317 return m->edges(*this);
318}
319
322{
323 return m->centre(*this, points);
324}
325
327inline Foam::scalar Foam::cellShape::mag(const UList<point>& points) const
328{
329 return m->mag(*this, points);
330}
331
332
333inline Foam::label Foam::cellShape::min() const
334{
335 label result = (empty() ? -1 : (*this)[0]);
336
337 for (const label pointi : *this)
338 {
339 if (result > pointi)
340 {
341 result = pointi;
343 }
344
345 return result;
346}
347
348
349inline Foam::label Foam::cellShape::max() const
350{
351 label result = (empty() ? -1 : (*this)[0]);
352
353 for (const label pointi : *this)
354 {
355 if (result < pointi)
356 {
357 result = pointi;
359 }
360
361 return result;
362}
363
364
365inline void Foam::cellShape::reset
366(
367 const cellModel& model,
368 const labelUList& labels,
369 const bool doCollapse
370)
371{
372 static_cast<labelList&>(*this) = labels;
373 m = &model;
374
375 if (doCollapse)
377 collapse();
378 }
379}
380
381
382template<unsigned N>
383inline void Foam::cellShape::reset
384(
385 const cellModel& model,
386 const FixedList<label, N>& labels,
387 const bool doCollapse
388)
389{
390 static_cast<labelList&>(*this) = labels;
391 m = &model;
392
393 if (doCollapse)
394 {
396 }
397}
398
399
400// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
401
402inline void Foam::cellShape::operator+=(const label vertexOffset)
403{
404 if (vertexOffset)
405 {
406 for (auto& vrt : static_cast<labelList&>(*this))
407 {
408 vrt += vertexOffset;
409 }
410 }
411}
412
413
414// ************************************************************************* //
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
A List with indirect addressing. Like IndirectList but does not store addressing.
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
bool empty() const noexcept
Definition UList.H:701
label size() const noexcept
Definition UList.H:706
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const Vector< Cmpt > & centre(const Foam::UList< Vector< Cmpt > > &) const noexcept
Return this (for point which is a typedef to Vector<scalar>).
Definition VectorI.H:67
scalar mag() const
The length (L2-norm) of the vector.
Definition VectorI.H:88
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
Maps a geometry to a set of cell primitives.
Definition cellModel.H:73
modelType
Enumeration of commonly used cellModel types.
Definition cellModel.H:82
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
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
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
label nFaces() const
Number of faces.
Definition cellShapeI.H:171
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
edge()
Default construct, with invalid vertex labels (-1).
Definition edgeI.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
constexpr face() noexcept=default
Default construct.
A class for handling words, derived from Foam::string.
Definition word.H:66
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
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
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299