Loading...
Searching...
No Matches
polyTopoChangeTemplates.C
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) 2017-2022 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 "polyMesh.H"
30#include "HashOps.H"
31#include "emptyPolyPatch.H"
32#include "mapPolyMesh.H"
33
34// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35
36template<class Type>
37void Foam::polyTopoChange::reorder
38(
39 const labelUList& oldToNew,
40 DynamicList<Type>& lst
41)
42{
43 // Create copy
44 DynamicList<Type> oldLst(lst);
45
46 forAll(oldLst, i)
47 {
48 const label newIdx = oldToNew[i];
49
50 if (newIdx >= 0)
51 {
52 lst[newIdx] = std::move(oldLst[i]);
53 }
54 }
55}
56
57
58template<class Type>
59void Foam::polyTopoChange::renumberKey
60(
61 const labelUList& oldToNew,
62 Map<Type>& map
63)
64{
65 Map<Type> newMap(map.capacity());
66
67 forAllConstIters(map, iter)
68 {
69 const label newKey = oldToNew[iter.key()];
70
71 if (newKey >= 0)
72 {
73 newMap.insert(newKey, iter.val());
74 }
75 }
76
77 map.transfer(newMap);
78}
79
80
81template<class Type>
83(
84 autoPtr<Type>& newMeshPtr,
85 const IOobject& io,
86 const polyMesh& mesh,
87 const labelUList& patchMap,
88 const bool syncParallel,
89 const bool orderCells,
90 const bool orderPoints
91)
92{
93 if (debug)
94 {
95 Pout<< "polyTopoChange::changeMesh"
96 << "(autoPtr<fvMesh>&, const IOobject&, const fvMesh&"
97 << ", const bool, const bool, const bool)"
98 << endl;
99 }
100
101 if (debug)
102 {
103 Pout<< "Old mesh:" << nl;
104 writeMeshStats(mesh, Pout);
105 }
106
107 // new mesh points
108 pointField newPoints;
109 // number of internal points
110 label nInternalPoints;
111 // patch slicing
112 labelList patchSizes;
113 labelList patchStarts;
114 // inflate maps
115 List<objectMap> pointsFromPoints;
116 List<objectMap> facesFromPoints;
117 List<objectMap> facesFromEdges;
118 List<objectMap> facesFromFaces;
119 List<objectMap> cellsFromPoints;
120 List<objectMap> cellsFromEdges;
121 List<objectMap> cellsFromFaces;
122 List<objectMap> cellsFromCells;
123
124 // old mesh info
125 List<Map<label>> oldPatchMeshPointMaps;
126 labelList oldPatchNMeshPoints;
127 labelList oldPatchStarts;
128 List<Map<label>> oldFaceZoneMeshPointMaps;
129
130 // Compact, reorder patch faces and calculate mesh/patch maps.
131 compactAndReorder
132 (
133 mesh,
134 patchMap, // from new to old patch
135 syncParallel,
136 orderCells,
137 orderPoints,
138
139 nInternalPoints,
140 newPoints,
141 patchSizes,
142 patchStarts,
143 pointsFromPoints,
144 facesFromPoints,
145 facesFromEdges,
146 facesFromFaces,
147 cellsFromPoints,
148 cellsFromEdges,
149 cellsFromFaces,
150 cellsFromCells,
151 oldPatchMeshPointMaps,
152 oldPatchNMeshPoints,
153 oldPatchStarts,
154 oldFaceZoneMeshPointMaps
155 );
156
157 const label nOldPoints(mesh.nPoints());
158 const label nOldFaces(mesh.nFaces());
159 const label nOldCells(mesh.nCells());
160 autoPtr<scalarField> oldCellVolumes(new scalarField(mesh.cellVolumes()));
161
162
163 // Create the mesh
164 // ~~~~~~~~~~~~~~~
165
166 //IOobject noReadIO(io);
167 //noReadIO.readOpt(IOobject::NO_READ);
168 //noReadIO.writeOpt(IOobject::AUTO_WRITE);
169 newMeshPtr.reset
170 (
171 new Type
172 (
173 io, //noReadIO
174 std::move(newPoints),
175 std::move(faces_),
176 std::move(faceOwner_),
177 std::move(faceNeighbour_)
178 )
179 );
180 Type& newMesh = *newMeshPtr;
181
182 // Clear out primitives
183 {
184 retiredPoints_.clearStorage();
185 region_.clearStorage();
186 }
187
188
189 if (debug)
190 {
191 // Some stats on changes
192 label nAdd, nInflate, nMerge, nRemove;
193 countMap(pointMap_, reversePointMap_, nAdd, nInflate, nMerge, nRemove);
194 Pout<< "Points:"
195 << " added(from point):" << nAdd
196 << " added(from nothing):" << nInflate
197 << " merged(into other point):" << nMerge
198 << " removed:" << nRemove
199 << nl;
200
201 countMap(faceMap_, reverseFaceMap_, nAdd, nInflate, nMerge, nRemove);
202 Pout<< "Faces:"
203 << " added(from face):" << nAdd
204 << " added(inflated):" << nInflate
205 << " merged(into other face):" << nMerge
206 << " removed:" << nRemove
207 << nl;
208
209 countMap(cellMap_, reverseCellMap_, nAdd, nInflate, nMerge, nRemove);
210 Pout<< "Cells:"
211 << " added(from cell):" << nAdd
212 << " added(inflated):" << nInflate
213 << " merged(into other cell):" << nMerge
214 << " removed:" << nRemove
215 << nl
216 << endl;
217 }
218
219
220 {
221 const polyBoundaryMesh& oldPatches = mesh.boundaryMesh();
222
223 polyPatchList newBoundary(patchMap.size());
224
225 forAll(patchMap, patchi)
226 {
227 const label oldPatchi = patchMap[patchi];
228
229 if (oldPatchi != -1)
230 {
231 newBoundary.set
232 (
233 patchi,
234 oldPatches[oldPatchi].clone
235 (
236 newMesh.boundaryMesh(),
237 patchi,
238 patchSizes[patchi],
239 patchStarts[patchi]
240 )
241 );
242 }
243 else
244 {
245 // Added patch
246 newBoundary.set
247 (
248 patchi,
250 (
251 "patch" + Foam::name(patchi),
252 patchSizes[patchi],
253 patchStarts[patchi],
254 patchi,
255 newMesh.boundaryMesh(),
257 )
258 );
259 }
260 }
261 newMesh.addFvPatches(newBoundary);
262 }
263
264
265 // Zones
266 // ~~~~~
267
268 // Start off from empty zones.
269 const pointZoneMesh& oldPointZones = mesh.pointZones();
270 List<pointZone*> pZonePtrs(oldPointZones.size());
271 {
272 forAll(oldPointZones, i)
273 {
274 pZonePtrs[i] = new pointZone
275 (
276 oldPointZones[i].name(),
277 i,
278 newMesh.pointZones()
279 );
280 }
281 }
282
283 const faceZoneMesh& oldFaceZones = mesh.faceZones();
284 List<faceZone*> fZonePtrs(oldFaceZones.size());
285 {
286 forAll(oldFaceZones, i)
287 {
288 fZonePtrs[i] = new faceZone
289 (
290 oldFaceZones[i].name(),
291 i,
292 newMesh.faceZones()
293 );
294 }
295 }
296
297 const cellZoneMesh& oldCellZones = mesh.cellZones();
298 List<cellZone*> cZonePtrs(oldCellZones.size());
299 {
300 forAll(oldCellZones, i)
301 {
302 cZonePtrs[i] = new cellZone
303 (
304 oldCellZones[i].name(),
305 i,
306 newMesh.cellZones()
307 );
308 }
309 }
310
311 newMesh.addZones(pZonePtrs, fZonePtrs, cZonePtrs);
312
313 // Inverse of point/face/cell zone addressing.
314 // For every preserved point/face/cells in zone give the old position.
315 // For added points, the index is set to -1
316 labelListList pointZoneMap(mesh.pointZones().size());
317 labelListList faceZoneFaceMap(mesh.faceZones().size());
318 labelListList cellZoneMap(mesh.cellZones().size());
319
320 resetZones(mesh, newMesh, pointZoneMap, faceZoneFaceMap, cellZoneMap);
321
322 // Clear zone info
323 {
324 pointZone_.clearStorage();
325 faceZone_.clearStorage();
326 faceZoneFlip_.clearStorage();
327 cellZone_.clearStorage();
328 }
329
330 // Patch point renumbering
331 // For every preserved point on a patch give the old position.
332 // For added points, the index is set to -1
333 labelListList patchPointMap(newMesh.boundaryMesh().size());
334 calcPatchPointMap
335 (
336 oldPatchMeshPointMaps,
337 patchMap,
338 newMesh.boundaryMesh(),
339 patchPointMap
340 );
341
342 // Create the face zone mesh point renumbering
343 labelListList faceZonePointMap(newMesh.faceZones().size());
344 calcFaceZonePointMap(newMesh, oldFaceZoneMeshPointMaps, faceZonePointMap);
345
346 if (debug)
347 {
348 Pout<< "New mesh:" << nl;
349 writeMeshStats(newMesh, Pout);
350 }
351
352 labelHashSet flipFaceFluxSet(HashSetOps::used(flipFaceFlux_));
353
355 (
356 newMesh,
357 nOldPoints,
358 nOldFaces,
359 nOldCells,
360
361 pointMap_,
362 pointsFromPoints,
363
364 faceMap_,
365 facesFromPoints,
366 facesFromEdges,
367 facesFromFaces,
368
369 cellMap_,
370 cellsFromPoints,
371 cellsFromEdges,
372 cellsFromFaces,
373 cellsFromCells,
374
375 reversePointMap_,
376 reverseFaceMap_,
377 reverseCellMap_,
378
379 flipFaceFluxSet,
380
381 patchPointMap,
382
383 pointZoneMap,
384
385 faceZonePointMap,
386 faceZoneFaceMap,
387 cellZoneMap,
388
389 newPoints, // if empty signals no inflation.
390 oldPatchStarts,
391 oldPatchNMeshPoints,
392 oldCellVolumes,
393 true // steal storage.
394 );
396 // At this point all member DynamicList (pointMap_, cellMap_ etc.) will
397 // be invalid.
398}
399
400
401template<class Type>
403(
404 autoPtr<Type>& newMeshPtr,
405 const IOobject& io,
406 const polyMesh& mesh,
407 const bool syncParallel,
408 const bool orderCells,
409 const bool orderPoints
410)
411{
412 return makeMesh
413 (
414 newMeshPtr,
415 io,
416 mesh,
417 identity(mesh.boundaryMesh().size()),
418 syncParallel,
419 orderCells,
420 orderPoints
421 );
422}
423
424
425// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A subset of mesh cells.
Definition cellZone.H:61
Empty front and back plane patch. Used for 2-D geometries.
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
A subset of mesh points.
Definition pointZone.H:64
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition polyMesh.H:663
autoPtr< mapPolyMesh > makeMesh(autoPtr< Type > &newMesh, const IOobject &io, const polyMesh &mesh, const labelUList &patchMap, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Create new mesh with old mesh patches. Additional dictionaries.
const scalarField & cellVolumes() const
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
static const word null
An empty word.
Definition word.H:84
dynamicFvMesh & mesh
const auto & io
auto & name
labelHashSet used(const bitSet &select)
Convert a bitset to a labelHashSet of the indices used.
Definition HashOps.C:26
Namespace for handling debugging switches.
Definition debug.C:45
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235