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-2021 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  List<polyPatch*> newBoundary(patchMap.size());
246 
247  forAll(patchMap, patchi)
248  {
249  const label oldPatchi = patchMap[patchi];
250 
251  if (oldPatchi != -1)
252  {
253  newBoundary[patchi] = oldPatches[oldPatchi].clone
254  (
255  newMesh.boundaryMesh(),
256  patchi,
257  patchSizes[patchi],
258  patchStarts[patchi]
259  ).ptr();
260  }
261  else
262  {
263  // Added patch
264  newBoundary[patchi] = new emptyPolyPatch
265  (
266  "patch" + Foam::name(patchi),
267  patchSizes[patchi],
268  patchStarts[patchi],
269  patchi,
270  newMesh.boundaryMesh(),
271  word::null
272  );
273  }
274  }
275  newMesh.addFvPatches(newBoundary);
276  }
277 
278 
279  // Zones
280  // ~~~~~
281 
282  // Start off from empty zones.
283  const pointZoneMesh& oldPointZones = mesh.pointZones();
284  List<pointZone*> pZonePtrs(oldPointZones.size());
285  {
286  forAll(oldPointZones, i)
287  {
288  pZonePtrs[i] = new pointZone
289  (
290  oldPointZones[i].name(),
291  i,
292  newMesh.pointZones()
293  );
294  }
295  }
296 
297  const faceZoneMesh& oldFaceZones = mesh.faceZones();
298  List<faceZone*> fZonePtrs(oldFaceZones.size());
299  {
300  forAll(oldFaceZones, i)
301  {
302  fZonePtrs[i] = new faceZone
303  (
304  oldFaceZones[i].name(),
305  i,
306  newMesh.faceZones()
307  );
308  }
309  }
310 
311  const cellZoneMesh& oldCellZones = mesh.cellZones();
312  List<cellZone*> cZonePtrs(oldCellZones.size());
313  {
314  forAll(oldCellZones, i)
315  {
316  cZonePtrs[i] = new cellZone
317  (
318  oldCellZones[i].name(),
319  i,
320  newMesh.cellZones()
321  );
322  }
323  }
324 
325  newMesh.addZones(pZonePtrs, fZonePtrs, cZonePtrs);
326 
327  // Inverse of point/face/cell zone addressing.
328  // For every preserved point/face/cells in zone give the old position.
329  // For added points, the index is set to -1
330  labelListList pointZoneMap(mesh.pointZones().size());
331  labelListList faceZoneFaceMap(mesh.faceZones().size());
332  labelListList cellZoneMap(mesh.cellZones().size());
333 
334  resetZones(mesh, newMesh, pointZoneMap, faceZoneFaceMap, cellZoneMap);
335 
336  // Clear zone info
337  {
338  pointZone_.clearStorage();
339  faceZone_.clearStorage();
340  faceZoneFlip_.clearStorage();
341  cellZone_.clearStorage();
342  }
343 
344  // Patch point renumbering
345  // For every preserved point on a patch give the old position.
346  // For added points, the index is set to -1
347  labelListList patchPointMap(newMesh.boundaryMesh().size());
348  calcPatchPointMap
349  (
350  oldPatchMeshPointMaps,
351  patchMap,
352  newMesh.boundaryMesh(),
353  patchPointMap
354  );
355 
356  // Create the face zone mesh point renumbering
357  labelListList faceZonePointMap(newMesh.faceZones().size());
358  calcFaceZonePointMap(newMesh, oldFaceZoneMeshPointMaps, faceZonePointMap);
359 
360  if (debug)
361  {
362  Pout<< "New mesh:" << nl;
363  writeMeshStats(newMesh, Pout);
364  }
365 
366  labelHashSet flipFaceFluxSet(HashSetOps::used(flipFaceFlux_));
367 
369  (
370  newMesh,
371  nOldPoints,
372  nOldFaces,
373  nOldCells,
374 
375  pointMap_,
376  pointsFromPoints,
377 
378  faceMap_,
379  facesFromPoints,
380  facesFromEdges,
381  facesFromFaces,
382 
383  cellMap_,
384  cellsFromPoints,
385  cellsFromEdges,
386  cellsFromFaces,
387  cellsFromCells,
388 
389  reversePointMap_,
390  reverseFaceMap_,
391  reverseCellMap_,
392 
393  flipFaceFluxSet,
394 
395  patchPointMap,
396 
397  pointZoneMap,
398 
399  faceZonePointMap,
400  faceZoneFaceMap,
401  cellZoneMap,
402 
403  newPoints, // if empty signals no inflation.
404  oldPatchStarts,
405  oldPatchNMeshPoints,
406  oldCellVolumes,
407  true // steal storage.
408  );
409 
410  // At this point all member DynamicList (pointMap_, cellMap_ etc.) will
411  // be invalid.
412 }
413 
414 
415 template<class Type>
417 (
418  autoPtr<Type>& newMeshPtr,
419  const IOobject& io,
420  const polyMesh& mesh,
421  const bool syncParallel,
422  const bool orderCells,
423  const bool orderPoints
424 )
425 {
426  return makeMesh
427  (
428  newMeshPtr,
429  io,
430  mesh,
431  identity(mesh.boundaryMesh().size()),
432  syncParallel,
433  orderCells,
434  orderPoints
435  );
436 }
437 
438 
439 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
mapPolyMesh.H
Foam::pointZone
A subset of mesh points.
Definition: pointZone.H:65
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::primitiveMesh::nPoints
label nPoints() const noexcept
Number of mesh points.
Definition: primitiveMeshI.H:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:62
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Field< vector >
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::ZoneMesh< pointZone, polyMesh >
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:486
Foam::PtrList::clone
PtrList< T > clone(Args &&... args) const
Make a copy by cloning each of the list elements.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::primitiveMesh::cellVolumes
const scalarField & cellVolumes() const
Definition: primitiveMeshCellCentresAndVols.C:96
emptyPolyPatch.H
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Foam::emptyPolyPatch
Empty front and back plane patch. Used for 2-D geometries.
Definition: emptyPolyPatch.H:50
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::List< label >
Foam::HashSetOps::used
labelHashSet used(const bitSet &select)
Convert a bitset to a labelHashSet of the indices used.
Definition: HashOps.C:35
Foam::UList< label >
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::polyTopoChange::makeMesh
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.
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:480
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
HashOps.H