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