All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
renumberMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  renumberMesh
26 
27 Description
28  Renumbers the cell list in order to reduce the bandwidth, reading and
29  renumbering all fields from all the time directories.
30 
31  By default uses bandCompression (CuthillMcKee) but will
32  read system/renumberMeshDict if -dict option is present
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "argList.H"
37 #include "IOobjectList.H"
38 #include "fvMesh.H"
39 #include "polyTopoChange.H"
40 #include "ReadFields.H"
41 #include "volFields.H"
42 #include "surfaceFields.H"
43 #include "SortableList.H"
44 #include "decompositionMethod.H"
45 #include "renumberMethod.H"
47 #include "CuthillMcKeeRenumber.H"
48 #include "fvMeshSubset.H"
49 #include "cellSet.H"
50 #include "faceSet.H"
51 #include "pointSet.H"
52 #include "systemDict.H"
53 
54 #ifdef FOAM_USE_ZOLTAN
55  #include "zoltanRenumber.H"
56 #endif
57 
58 
59 using namespace Foam;
60 
61 
62 // Create named field from labelList for postprocessing
63 tmp<volScalarField> createScalarField
64 (
65  const fvMesh& mesh,
66  const word& name,
67  const labelList& elems
68 )
69 {
71  (
73  (
74  name,
75  mesh,
77  zeroGradientFvPatchScalarField::typeName
78  )
79  );
80  volScalarField& fld = tfld.ref();
81 
82  forAll(fld, celli)
83  {
84  fld[celli] = elems[celli];
85  }
86 
87  return tfld;
88 }
89 
90 
91 // Calculate band of matrix
92 label getBand(const labelList& owner, const labelList& neighbour)
93 {
94  label band = 0;
95 
96  forAll(neighbour, facei)
97  {
98  label diff = neighbour[facei] - owner[facei];
99 
100  if (diff > band)
101  {
102  band = diff;
103  }
104  }
105  return band;
106 }
107 
108 
109 // Calculate band of matrix
110 void getBand
111 (
112  const bool calculateIntersect,
113  const label nCells,
114  const labelList& owner,
115  const labelList& neighbour,
116  label& bandwidth,
117  scalar& profile, // scalar to avoid overflow
118  scalar& sumSqrIntersect // scalar to avoid overflow
119 )
120 {
121  labelList cellBandwidth(nCells, 0);
122  scalarField nIntersect(nCells, 0.0);
123 
124  forAll(neighbour, facei)
125  {
126  label own = owner[facei];
127  label nei = neighbour[facei];
128 
129  // Note: mag not necessary for correct (upper-triangular) ordering.
130  label diff = nei-own;
131  cellBandwidth[nei] = max(cellBandwidth[nei], diff);
132  }
133 
134  bandwidth = max(cellBandwidth);
135 
136  // Do not use field algebra because of conversion label to scalar
137  profile = 0.0;
138  forAll(cellBandwidth, celli)
139  {
140  profile += 1.0*cellBandwidth[celli];
141  }
142 
143  sumSqrIntersect = 0.0;
144  if (calculateIntersect)
145  {
146  forAll(nIntersect, celli)
147  {
148  for (label colI = celli-cellBandwidth[celli]; colI <= celli; colI++)
149  {
150  nIntersect[colI] += 1.0;
151  }
152  }
153 
154  sumSqrIntersect = sum(Foam::sqr(nIntersect));
155  }
156 }
157 
158 
159 // Determine upper-triangular face order
160 labelList getFaceOrder
161 (
162  const primitiveMesh& mesh,
163  const labelList& cellOrder // New to old cell
164 )
165 {
166  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
167 
168  labelList oldToNewFace(mesh.nFaces(), -1);
169 
170  label newFacei = 0;
171 
172  labelList nbr;
173  labelList order;
174 
175  forAll(cellOrder, newCelli)
176  {
177  label oldCelli = cellOrder[newCelli];
178 
179  const cell& cFaces = mesh.cells()[oldCelli];
180 
181  // Neighbouring cells
182  nbr.setSize(cFaces.size());
183 
184  forAll(cFaces, i)
185  {
186  label facei = cFaces[i];
187 
188  if (mesh.isInternalFace(facei))
189  {
190  // Internal face. Get cell on other side.
191  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
192  if (nbrCelli == newCelli)
193  {
194  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
195  }
196 
197  if (newCelli < nbrCelli)
198  {
199  // Celli is master
200  nbr[i] = nbrCelli;
201  }
202  else
203  {
204  // nbrCell is master. Let it handle this face.
205  nbr[i] = -1;
206  }
207  }
208  else
209  {
210  // External face. Do later.
211  nbr[i] = -1;
212  }
213  }
214 
215  order.setSize(nbr.size());
216  sortedOrder(nbr, order);
217 
218  forAll(order, i)
219  {
220  label index = order[i];
221  if (nbr[index] != -1)
222  {
223  oldToNewFace[cFaces[index]] = newFacei++;
224  }
225  }
226  }
227 
228  // Leave patch faces intact.
229  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
230  {
231  oldToNewFace[facei] = facei;
232  }
233 
234 
235  // Check done all faces.
236  forAll(oldToNewFace, facei)
237  {
238  if (oldToNewFace[facei] == -1)
239  {
241  << "Did not determine new position" << " for face " << facei
242  << abort(FatalError);
243  }
244  }
245 
246  return invert(mesh.nFaces(), oldToNewFace);
247 }
248 
249 
250 // Determine face order such that inside region faces are sorted
251 // upper-triangular but in between region faces are handled like boundary faces.
252 labelList getRegionFaceOrder
253 (
254  const primitiveMesh& mesh,
255  const labelList& cellOrder, // New to old cell
256  const labelList& cellToRegion // Old cell to region
257 )
258 {
259  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
260 
261  labelList oldToNewFace(mesh.nFaces(), -1);
262 
263  label newFacei = 0;
264 
265  label prevRegion = -1;
266 
267  forAll(cellOrder, newCelli)
268  {
269  label oldCelli = cellOrder[newCelli];
270 
271  if (cellToRegion[oldCelli] != prevRegion)
272  {
273  prevRegion = cellToRegion[oldCelli];
274  }
275 
276  const cell& cFaces = mesh.cells()[oldCelli];
277 
278  SortableList<label> nbr(cFaces.size());
279 
280  forAll(cFaces, i)
281  {
282  label facei = cFaces[i];
283 
284  if (mesh.isInternalFace(facei))
285  {
286  // Internal face. Get cell on other side.
287  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
288  if (nbrCelli == newCelli)
289  {
290  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
291  }
292 
293  if (cellToRegion[oldCelli] != cellToRegion[cellOrder[nbrCelli]])
294  {
295  // Treat like external face. Do later.
296  nbr[i] = -1;
297  }
298  else if (newCelli < nbrCelli)
299  {
300  // Celli is master
301  nbr[i] = nbrCelli;
302  }
303  else
304  {
305  // nbrCell is master. Let it handle this face.
306  nbr[i] = -1;
307  }
308  }
309  else
310  {
311  // External face. Do later.
312  nbr[i] = -1;
313  }
314  }
315 
316  nbr.sort();
317 
318  forAll(nbr, i)
319  {
320  if (nbr[i] != -1)
321  {
322  oldToNewFace[cFaces[nbr.indices()[i]]] = newFacei++;
323  }
324  }
325  }
326 
327  // Do region interfaces
328  label nRegions = max(cellToRegion)+1;
329  {
330  // Sort in increasing region
331  SortableList<label> sortKey(mesh.nFaces(), labelMax);
332 
333  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
334  {
335  label ownRegion = cellToRegion[mesh.faceOwner()[facei]];
336  label neiRegion = cellToRegion[mesh.faceNeighbour()[facei]];
337 
338  if (ownRegion != neiRegion)
339  {
340  sortKey[facei] =
341  min(ownRegion, neiRegion)*nRegions
342  +max(ownRegion, neiRegion);
343  }
344  }
345  sortKey.sort();
346 
347  // Extract.
348  label prevKey = -1;
349  forAll(sortKey, i)
350  {
351  label key = sortKey[i];
352 
353  if (key == labelMax)
354  {
355  break;
356  }
357 
358  if (prevKey != key)
359  {
360  prevKey = key;
361  }
362 
363  oldToNewFace[sortKey.indices()[i]] = newFacei++;
364  }
365  }
366 
367  // Leave patch faces intact.
368  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
369  {
370  oldToNewFace[facei] = facei;
371  }
372 
373 
374  // Check done all faces.
375  forAll(oldToNewFace, facei)
376  {
377  if (oldToNewFace[facei] == -1)
378  {
380  << "Did not determine new position"
381  << " for face " << facei
382  << abort(FatalError);
383  }
384  }
385 
386  return invert(mesh.nFaces(), oldToNewFace);
387 }
388 
389 
390 // cellOrder: old cell for every new cell
391 // faceOrder: old face for every new face. Ordering of boundary faces not
392 // changed.
393 autoPtr<polyTopoChangeMap> reorderMesh
394 (
395  polyMesh& mesh,
396  const labelList& cellOrder,
397  const labelList& faceOrder
398 )
399 {
400  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
401  labelList reverseFaceOrder(invert(faceOrder.size(), faceOrder));
402 
403  faceList newFaces(reorder(reverseFaceOrder, mesh.faces()));
404  labelList newOwner
405  (
406  renumber
407  (
408  reverseCellOrder,
409  reorder(reverseFaceOrder, mesh.faceOwner())
410  )
411  );
412  labelList newNeighbour
413  (
414  renumber
415  (
416  reverseCellOrder,
417  reorder(reverseFaceOrder, mesh.faceNeighbour())
418  )
419  );
420 
421  // Check if any faces need swapping.
422  labelHashSet flipFaceFlux(newOwner.size());
423  forAll(newNeighbour, facei)
424  {
425  label own = newOwner[facei];
426  label nei = newNeighbour[facei];
427 
428  if (nei < own)
429  {
430  newFaces[facei].flip();
431  Swap(newOwner[facei], newNeighbour[facei]);
432  flipFaceFlux.insert(facei);
433  }
434  }
435 
436  const polyBoundaryMesh& patches = mesh.boundaryMesh();
437  labelList patchSizes(patches.size());
438  labelList patchStarts(patches.size());
439  labelList oldPatchNMeshPoints(patches.size());
440  labelListList patchPointMap(patches.size());
441 
442  forAll(patches, patchi)
443  {
444  patchSizes[patchi] = patches[patchi].size();
445  patchStarts[patchi] = patches[patchi].start();
446  oldPatchNMeshPoints[patchi] = patches[patchi].nPoints();
447  patchPointMap[patchi] = identity(patches[patchi].nPoints());
448  }
449 
450  mesh.resetPrimitives
451  (
452  NullObjectMove<pointField>(),
453  move(newFaces),
454  move(newOwner),
455  move(newNeighbour),
456  patchSizes,
457  patchStarts,
458  true
459  );
460 
461 
462  // Re-do the faceZones
463  {
464  meshFaceZones& faceZones = mesh.faceZones();
465  faceZones.clearAddressing();
466  forAll(faceZones, zoneI)
467  {
468  faceZone& fZone = faceZones[zoneI];
469  labelList newAddressing(fZone.size());
470  boolList newFlipMap(fZone.size());
471  forAll(fZone, i)
472  {
473  label oldFacei = fZone[i];
474  newAddressing[i] = reverseFaceOrder[oldFacei];
475  if (flipFaceFlux.found(newAddressing[i]))
476  {
477  newFlipMap[i] = !fZone.flipMap()[i];
478  }
479  else
480  {
481  newFlipMap[i] = fZone.flipMap()[i];
482  }
483  }
484  labelList newToOld;
485  sortedOrder(newAddressing, newToOld);
486  fZone.resetAddressing
487  (
488  UIndirectList<label>(newAddressing, newToOld)(),
489  UIndirectList<bool>(newFlipMap, newToOld)()
490  );
491  }
492  }
493  // Re-do the cellZones
494  {
495  meshCellZones& cellZones = mesh.cellZones();
496  cellZones.clearAddressing();
497  forAll(cellZones, zoneI)
498  {
499  cellZones[zoneI] = UIndirectList<label>
500  (
501  reverseCellOrder,
502  cellZones[zoneI]
503  )();
504  Foam::sort(cellZones[zoneI]);
505  }
506  }
507 
508 
510  (
512  (
513  mesh, // const polyMesh& mesh,
514  mesh.nPoints(), // nOldPoints,
515  mesh.nFaces(), // nOldFaces,
516  mesh.nCells(), // nOldCells,
517  identity(mesh.nPoints()), // pointMap,
518  List<objectMap>(0), // pointsFromPoints,
519  faceOrder, // faceMap,
520  List<objectMap>(0), // facesFromPoints,
521  List<objectMap>(0), // facesFromEdges,
522  List<objectMap>(0), // facesFromFaces,
523  cellOrder, // cellMap,
524  List<objectMap>(0), // cellsFromPoints,
525  List<objectMap>(0), // cellsFromEdges,
526  List<objectMap>(0), // cellsFromFaces,
527  List<objectMap>(0), // cellsFromCells,
528  identity(mesh.nPoints()), // reversePointMap,
529  reverseFaceOrder, // reverseFaceMap,
530  reverseCellOrder, // reverseCellMap,
531  flipFaceFlux, // flipFaceFlux,
532  patchPointMap, // patchPointMap,
533  labelListList(0), // pointZoneMap,
534  labelListList(0), // faceZonePointMap,
535  labelListList(0), // faceZoneFaceMap,
536  labelListList(0), // cellZoneMap,
537  pointField(0), // preMotionPoints,
538  patchStarts, // oldPatchStarts,
539  oldPatchNMeshPoints, // oldPatchNMeshPoints
540  autoPtr<scalarField>() // oldCellVolumes
541  )
542  );
543 }
544 
545 
546 // Return new to old cell numbering
547 labelList regionRenumber
548 (
549  const renumberMethod& method,
550  const fvMesh& mesh,
551  const labelList& cellToRegion
552 )
553 {
554  Info<< "Determining cell order:" << endl;
555 
556  labelList cellOrder(cellToRegion.size());
557 
558  label nRegions = max(cellToRegion)+1;
559 
560  labelListList regionToCells(invertOneToMany(nRegions, cellToRegion));
561 
562  label celli = 0;
563 
564  forAll(regionToCells, regionI)
565  {
566  Info<< " region " << regionI << " starts at " << celli << endl;
567 
568  // Make sure no parallel comms
569  bool oldParRun = UPstream::parRun();
570  UPstream::parRun() = false;
571 
572  // Per region do a reordering.
573  fvMeshSubset subsetter(mesh);
574  subsetter.setLargeCellSubset(cellToRegion, regionI);
575 
576  const fvMesh& subMesh = subsetter.subMesh();
577 
578  labelList subCellOrder = method.renumber
579  (
580  subMesh,
581  subMesh.cellCentres()
582  );
583 
584  // Restore state
585  UPstream::parRun() = oldParRun;
586 
587  const labelList& cellMap = subsetter.cellMap();
588 
589  forAll(subCellOrder, i)
590  {
591  cellOrder[celli++] = cellMap[subCellOrder[i]];
592  }
593  }
594  Info<< endl;
595 
596  return cellOrder;
597 }
598 
599 
600 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
601 
602 int main(int argc, char *argv[])
603 {
605  (
606  "Renumber mesh to minimise bandwidth"
607  );
608 
609  #include "addRegionOption.H"
610  #include "addOverwriteOption.H"
611  #include "addTimeOptions.H"
612  #include "addDictOption.H"
614  (
615  "frontWidth",
616  "calculate the rms of the frontwidth"
617  );
619  (
620  "noFields",
621  "do not update fields"
622  );
623 
624  #include "setRootCase.H"
625  #include "createTime.H"
626  runTime.functionObjects().off();
627 
628  // Force linker to include zoltan symbols. This section is only needed since
629  // Zoltan is a static library
630  #ifdef FOAM_USE_ZOLTAN
631  Info<< "renumberMesh built with zoltan support." << nl << endl;
632  (void)zoltanRenumber::typeName;
633  #endif
634 
635  // Get times list
636  instantList Times = runTime.times();
637 
638  // Set startTime and endTime depending on -time and -latestTime options
639  #include "checkTimeOptions.H"
640 
641  runTime.setTime(Times[startTime], startTime);
642 
643  #include "createNamedMesh.H"
644  const word oldInstance = mesh.pointsInstance();
645 
646  const bool readDict = args.optionFound("dict");
647  const bool doFrontWidth = args.optionFound("frontWidth");
648  const bool overwrite = args.optionFound("overwrite");
649  const bool fields = !args.optionFound("noFields");
650 
651  label band;
652  scalar profile;
653  scalar sumSqrIntersect;
654  getBand
655  (
656  doFrontWidth,
657  mesh.nCells(),
658  mesh.faceOwner(),
659  mesh.faceNeighbour(),
660  band,
661  profile,
662  sumSqrIntersect
663  );
664 
665  reduce(band, maxOp<label>());
666  reduce(profile, sumOp<scalar>());
667  scalar rmsFrontwidth = Foam::sqrt
668  (
670  (
671  sumSqrIntersect,
672  sumOp<scalar>()
673  )/mesh.globalData().nTotalCells()
674  );
675 
676  Info<< "Mesh size: " << mesh.globalData().nTotalCells() << nl
677  << "Before renumbering :" << nl
678  << " band : " << band << nl
679  << " profile : " << profile << nl;
680 
681  if (doFrontWidth)
682  {
683  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
684  }
685 
686  Info<< endl;
687 
688  bool sortCoupledFaceCells = false;
689  bool writeMaps = false;
690  bool orderPoints = false;
691  label blockSize = 0;
692  bool renumberSets = true;
693 
694  // Construct renumberMethod
695  autoPtr<dictionary> renumberDictPtr;
696  autoPtr<renumberMethod> renumberPtr;
697 
698  if (readDict)
699  {
700  renumberDictPtr.reset
701  (
702  new dictionary(systemDict("renumberMeshDict", args, mesh))
703  );
704  const dictionary& renumberDict = renumberDictPtr();
705 
706  renumberPtr = renumberMethod::New(renumberDict);
707 
708  sortCoupledFaceCells = renumberDict.lookupOrDefault
709  (
710  "sortCoupledFaceCells",
711  false
712  );
713  if (sortCoupledFaceCells)
714  {
715  Info<< "Sorting cells on coupled boundaries to be last." << nl
716  << endl;
717  }
718 
719  blockSize = renumberDict.lookupOrDefault("blockSize", 0);
720  if (blockSize > 0)
721  {
722  Info<< "Ordering cells into regions of size " << blockSize
723  << " (using decomposition);"
724  << " ordering faces into region-internal and region-external."
725  << nl << endl;
726 
727  if (blockSize < 0 || blockSize >= mesh.nCells())
728  {
730  << "Block size " << blockSize
731  << " should be positive integer"
732  << " and less than the number of cells in the mesh."
733  << exit(FatalError);
734  }
735  }
736 
737  orderPoints = renumberDict.lookupOrDefault("orderPoints", false);
738  if (orderPoints)
739  {
740  Info<< "Ordering points into internal and boundary points." << nl
741  << endl;
742  }
743 
744  renumberDict.lookup("writeMaps") >> writeMaps;
745  if (writeMaps)
746  {
747  Info<< "Writing renumber maps (new to old) to polyMesh." << nl
748  << endl;
749  }
750 
751  renumberSets = renumberDict.lookupOrDefault("renumberSets", true);
752  }
753  else
754  {
755  Info<< "Using default renumberMethod." << nl << endl;
756  dictionary renumberDict;
757  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
758  }
759 
760  Info<< "Selecting renumberMethod " << renumberPtr().type() << nl << endl;
761 
762 
763 
764  // Read parallel reconstruct maps
765  labelIOList cellProcAddressing
766  (
767  IOobject
768  (
769  "cellProcAddressing",
770  mesh.facesInstance(),
772  mesh,
774  ),
775  labelList(0)
776  );
777 
778  labelIOList faceProcAddressing
779  (
780  IOobject
781  (
782  "faceProcAddressing",
783  mesh.facesInstance(),
785  mesh,
787  ),
788  labelList(0)
789  );
790  labelIOList pointProcAddressing
791  (
792  IOobject
793  (
794  "pointProcAddressing",
795  mesh.pointsInstance(),
797  mesh,
799  ),
800  labelList(0)
801  );
802 
803 
804  // Read objects in time directory
805  IOobjectList objects(mesh, runTime.timeName());
806 
807  if (fields) Info<< "Reading geometric fields" << nl << endl;
808 
809  #include "readVolFields.H"
810  #include "readSurfaceFields.H"
811  #include "readPointFields.H"
812 
813  Info<< endl;
814 
815  // From renumbering:
816  // - from new cell/face back to original cell/face
817  labelList cellOrder;
818  labelList faceOrder;
819 
820  if (blockSize > 0)
821  {
822  // Renumbering in two phases. Should be done in one so mapping of
823  // fields is done correctly!
824 
825  label nBlocks = mesh.nCells()/blockSize;
826  Info<< "nBlocks = " << nBlocks << endl;
827 
828  // Read decompositionMethod dictionary
829  dictionary decomposeDict(renumberDictPtr().subDict("blockCoeffs"));
830  decomposeDict.set("numberOfSubdomains", nBlocks);
831 
832  bool oldParRun = UPstream::parRun();
833  UPstream::parRun() = false;
834 
835  autoPtr<decompositionMethod> decomposePtr
836  (
838  );
839 
840  labelList cellToRegion
841  (
842  decomposePtr().decompose
843  (
844  mesh,
845  mesh.cellCentres()
846  )
847  );
848 
849  // Restore state
850  UPstream::parRun() = oldParRun;
851 
852  // For debugging: write out region
853  createScalarField
854  (
855  mesh,
856  "cellDist",
857  cellToRegion
858  )().write();
859 
860  Info<< nl << "Written decomposition as volScalarField to "
861  << "cellDist for use in postprocessing."
862  << nl << endl;
863 
864 
865  cellOrder = regionRenumber(renumberPtr(), mesh, cellToRegion);
866 
867  // Determine new to old face order with new cell numbering
868  faceOrder = getRegionFaceOrder
869  (
870  mesh,
871  cellOrder,
872  cellToRegion
873  );
874  }
875  else
876  {
877  // Determines sorted back to original cell ordering
878  cellOrder = renumberPtr().renumber
879  (
880  mesh,
881  mesh.cellCentres()
882  );
883 
884  if (sortCoupledFaceCells)
885  {
886  // Change order so all coupled patch faceCells are at the end.
887  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
888 
889  // Collect all boundary cells on coupled patches
890  label nBndCells = 0;
891  forAll(pbm, patchi)
892  {
893  if (pbm[patchi].coupled())
894  {
895  nBndCells += pbm[patchi].size();
896  }
897  }
898 
899  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
900 
901  labelList bndCells(nBndCells);
902  labelList bndCellMap(nBndCells);
903  nBndCells = 0;
904  forAll(pbm, patchi)
905  {
906  if (pbm[patchi].coupled())
907  {
908  const labelUList& faceCells = pbm[patchi].faceCells();
909  forAll(faceCells, i)
910  {
911  label celli = faceCells[i];
912 
913  if (reverseCellOrder[celli] != -1)
914  {
915  bndCells[nBndCells] = celli;
916  bndCellMap[nBndCells++] = reverseCellOrder[celli];
917  reverseCellOrder[celli] = -1;
918  }
919  }
920  }
921  }
922  bndCells.setSize(nBndCells);
923  bndCellMap.setSize(nBndCells);
924 
925  // Sort
926  labelList order;
927  sortedOrder(bndCellMap, order);
928 
929  // Redo newReverseCellOrder
930  labelList newReverseCellOrder(mesh.nCells(), -1);
931 
932  label sortedI = mesh.nCells();
933  forAllReverse(order, i)
934  {
935  label origCelli = bndCells[order[i]];
936  newReverseCellOrder[origCelli] = --sortedI;
937  }
938 
939  Info<< "Ordered all " << nBndCells << " cells with a coupled face"
940  << " to the end of the cell list, starting at " << sortedI
941  << endl;
942 
943  // Compact
944  sortedI = 0;
945  forAll(cellOrder, newCelli)
946  {
947  label origCelli = cellOrder[newCelli];
948  if (newReverseCellOrder[origCelli] == -1)
949  {
950  newReverseCellOrder[origCelli] = sortedI++;
951  }
952  }
953 
954  // Update sorted back to original (unsorted) map
955  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
956  }
957 
958 
959  // Determine new to old face order with new cell numbering
960  faceOrder = getFaceOrder
961  (
962  mesh,
963  cellOrder // New to old cell
964  );
965  }
966 
967 
968  if (!overwrite)
969  {
970  runTime++;
971  }
972 
973 
974  // Change the mesh.
975  autoPtr<polyTopoChangeMap> map = reorderMesh(mesh, cellOrder, faceOrder);
976 
977 
978  if (orderPoints)
979  {
980  polyTopoChange meshMod(mesh);
981  autoPtr<polyTopoChangeMap> pointOrderMap = meshMod.changeMesh
982  (
983  mesh,
984  false, // inflate
985  true, // syncParallel
986  false, // orderCells
987  orderPoints // orderPoints
988  );
989 
990  // Combine point reordering into map.
991  const_cast<labelList&>(map().pointMap()) = UIndirectList<label>
992  (
993  map().pointMap(),
994  pointOrderMap().pointMap()
995  )();
996 
998  (
999  pointOrderMap().reversePointMap(),
1000  const_cast<labelList&>(map().reversePointMap())
1001  );
1002  }
1003 
1004 
1005  // Update fields
1006  mesh.topoChange(map);
1007 
1008  // Update proc maps
1009  if
1010  (
1011  cellProcAddressing.headerOk()
1012  && cellProcAddressing.size() == mesh.nCells()
1013  )
1014  {
1015  Info<< "Renumbering processor cell decomposition map "
1016  << cellProcAddressing.name() << endl;
1017 
1018  cellProcAddressing = labelList
1019  (
1020  UIndirectList<label>(cellProcAddressing, map().cellMap())
1021  );
1022  }
1023  if
1024  (
1025  faceProcAddressing.headerOk()
1026  && faceProcAddressing.size() == mesh.nFaces()
1027  )
1028  {
1029  Info<< "Renumbering processor face decomposition map "
1030  << faceProcAddressing.name() << endl;
1031 
1032  faceProcAddressing = labelList
1033  (
1034  UIndirectList<label>(faceProcAddressing, map().faceMap())
1035  );
1036 
1037  // Detect any flips.
1038  const labelHashSet& fff = map().flipFaceFlux();
1039  forAllConstIter(labelHashSet, fff, iter)
1040  {
1041  label facei = iter.key();
1042  label masterFacei = faceProcAddressing[facei];
1043 
1044  faceProcAddressing[facei] = -masterFacei;
1045 
1046  if (masterFacei == 0)
1047  {
1049  << " masterFacei:" << masterFacei << exit(FatalError);
1050  }
1051  }
1052  }
1053  if
1054  (
1055  pointProcAddressing.headerOk()
1056  && pointProcAddressing.size() == mesh.nPoints()
1057  )
1058  {
1059  Info<< "Renumbering processor point decomposition map "
1060  << pointProcAddressing.name() << endl;
1061 
1062  pointProcAddressing = labelList
1063  (
1064  UIndirectList<label>(pointProcAddressing, map().pointMap())
1065  );
1066  }
1067 
1068 
1069  // Move mesh (since morphing might not do this)
1070  if (map().hasMotionPoints())
1071  {
1072  mesh.setPoints(map().preMotionPoints());
1073  }
1074 
1075 
1076  {
1077  label band;
1078  scalar profile;
1079  scalar sumSqrIntersect;
1080  getBand
1081  (
1082  doFrontWidth,
1083  mesh.nCells(),
1084  mesh.faceOwner(),
1085  mesh.faceNeighbour(),
1086  band,
1087  profile,
1088  sumSqrIntersect
1089  );
1090  reduce(band, maxOp<label>());
1091  reduce(profile, sumOp<scalar>());
1092  scalar rmsFrontwidth = Foam::sqrt
1093  (
1094  returnReduce
1095  (
1096  sumSqrIntersect,
1097  sumOp<scalar>()
1098  )/mesh.globalData().nTotalCells()
1099  );
1100 
1101  Info<< "After renumbering :" << nl
1102  << " band : " << band << nl
1103  << " profile : " << profile << nl;
1104 
1105  if (doFrontWidth)
1106  {
1107 
1108  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1109  }
1110 
1111  Info<< endl;
1112  }
1113 
1114  if (orderPoints)
1115  {
1116  // Force edge calculation (since only reason that points would need to
1117  // be sorted)
1118  (void)mesh.edges();
1119 
1120  label nTotPoints = returnReduce
1121  (
1122  mesh.nPoints(),
1123  sumOp<label>()
1124  );
1125  label nTotIntPoints = returnReduce
1126  (
1127  mesh.nInternalPoints(),
1128  sumOp<label>()
1129  );
1130 
1131  label nTotEdges = returnReduce
1132  (
1133  mesh.nEdges(),
1134  sumOp<label>()
1135  );
1136  label nTotIntEdges = returnReduce
1137  (
1138  mesh.nInternalEdges(),
1139  sumOp<label>()
1140  );
1141  label nTotInt0Edges = returnReduce
1142  (
1143  mesh.nInternal0Edges(),
1144  sumOp<label>()
1145  );
1146  label nTotInt1Edges = returnReduce
1147  (
1148  mesh.nInternal1Edges(),
1149  sumOp<label>()
1150  );
1151 
1152  Info<< "Points:" << nl
1153  << " total : " << nTotPoints << nl
1154  << " internal: " << nTotIntPoints << nl
1155  << " boundary: " << nTotPoints-nTotIntPoints << nl
1156  << "Edges:" << nl
1157  << " total : " << nTotEdges << nl
1158  << " internal: " << nTotIntEdges << nl
1159  << " internal using 0 boundary points: "
1160  << nTotInt0Edges << nl
1161  << " internal using 1 boundary points: "
1162  << nTotInt1Edges-nTotInt0Edges << nl
1163  << " internal using 2 boundary points: "
1164  << nTotIntEdges-nTotInt1Edges << nl
1165  << " boundary: " << nTotEdges-nTotIntEdges << nl
1166  << endl;
1167  }
1168 
1169 
1170  if (overwrite)
1171  {
1172  mesh.setInstance(oldInstance);
1173  }
1174 
1175  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1176 
1177  mesh.write();
1178 
1179  if (cellProcAddressing.headerOk())
1180  {
1181  cellProcAddressing.instance() = mesh.facesInstance();
1182  if (cellProcAddressing.size() == mesh.nCells())
1183  {
1184  cellProcAddressing.write();
1185  }
1186  else
1187  {
1188  // procAddressing file no longer valid. Delete it.
1189  const fileName fName(cellProcAddressing.filePath());
1190  if (fName.size())
1191  {
1192  Info<< "Deleting inconsistent processor cell decomposition"
1193  << " map " << fName << endl;
1194  rm(fName);
1195  }
1196  }
1197  }
1198 
1199  if (faceProcAddressing.headerOk())
1200  {
1201  faceProcAddressing.instance() = mesh.facesInstance();
1202  if (faceProcAddressing.size() == mesh.nFaces())
1203  {
1204  faceProcAddressing.write();
1205  }
1206  else
1207  {
1208  const fileName fName(faceProcAddressing.filePath());
1209  if (fName.size())
1210  {
1211  Info<< "Deleting inconsistent processor face decomposition"
1212  << " map " << fName << endl;
1213  rm(fName);
1214  }
1215  }
1216  }
1217 
1218  if (pointProcAddressing.headerOk())
1219  {
1220  pointProcAddressing.instance() = mesh.facesInstance();
1221  if (pointProcAddressing.size() == mesh.nPoints())
1222  {
1223  pointProcAddressing.write();
1224  }
1225  else
1226  {
1227  const fileName fName(pointProcAddressing.filePath());
1228  if (fName.size())
1229  {
1230  Info<< "Deleting inconsistent processor point decomposition"
1231  << " map " << fName << endl;
1232  rm(fName);
1233  }
1234  }
1235  }
1236 
1237  if (writeMaps)
1238  {
1239  // For debugging: write out region
1240  createScalarField
1241  (
1242  mesh,
1243  "origCellID",
1244  map().cellMap()
1245  )().write();
1246 
1247  createScalarField
1248  (
1249  mesh,
1250  "cellID",
1251  identity(mesh.nCells())
1252  )().write();
1253 
1254  Info<< nl << "Written current cellID and origCellID as volScalarField"
1255  << " for use in postprocessing."
1256  << nl << endl;
1257 
1258  labelIOList
1259  (
1260  IOobject
1261  (
1262  "cellMap",
1263  mesh.facesInstance(),
1265  mesh,
1268  false
1269  ),
1270  map().cellMap()
1271  ).write();
1272 
1273  labelIOList
1274  (
1275  IOobject
1276  (
1277  "faceMap",
1278  mesh.facesInstance(),
1280  mesh,
1283  false
1284  ),
1285  map().faceMap()
1286  ).write();
1287 
1288  labelIOList
1289  (
1290  IOobject
1291  (
1292  "pointMap",
1293  mesh.facesInstance(),
1295  mesh,
1298  false
1299  ),
1300  map().pointMap()
1301  ).write();
1302  }
1303 
1304 
1305  // Renumber sets if required
1306  if (renumberSets)
1307  {
1308  Info<< endl;
1309 
1310  // Read sets
1311  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1312 
1313  {
1314  IOobjectList cSets(objects.lookupClass(cellSet::typeName));
1315  if (cSets.size())
1316  {
1317  Info<< "Renumbering cellSets:" << endl;
1318  forAllConstIter(IOobjectList, cSets, iter)
1319  {
1320  cellSet cs(*iter());
1321  Info<< " " << cs.name() << endl;
1322  cs.topoChange(map());
1323  cs.instance() = mesh.facesInstance();
1324  cs.write();
1325  }
1326  }
1327  }
1328 
1329  {
1330  IOobjectList fSets(objects.lookupClass(faceSet::typeName));
1331  if (fSets.size())
1332  {
1333  Info<< "Renumbering faceSets:" << endl;
1334  forAllConstIter(IOobjectList, fSets, iter)
1335  {
1336  faceSet fs(*iter());
1337  Info<< " " << fs.name() << endl;
1338  fs.topoChange(map());
1339  fs.instance() = mesh.facesInstance();
1340  fs.write();
1341  }
1342  }
1343  }
1344 
1345  {
1346  IOobjectList pSets(objects.lookupClass(pointSet::typeName));
1347  if (pSets.size())
1348  {
1349  Info<< "Renumbering pointSets:" << endl;
1350  forAllConstIter(IOobjectList, pSets, iter)
1351  {
1352  pointSet ps(*iter());
1353  Info<< " " << ps.name() << endl;
1354  ps.topoChange(map());
1355  ps.instance() = mesh.facesInstance();
1356  ps.write();
1357  }
1358  }
1359  }
1360  }
1361 
1362  Info<< "End\n" << endl;
1363 
1364  return 0;
1365 }
1366 
1367 
1368 // ************************************************************************* //
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:403
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual Ostream & write(const char)=0
Write character.
virtual void resetAddressing(const labelUList &, const boolList &)
Reset addressing and flip map (clearing demand-driven data)
Definition: faceZone.C:379
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
A class for handling file names.
Definition: fileName.H:79
A list of face labels.
Definition: faceSet.H:48
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:888
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
A set of point labels.
Definition: pointSet.H:48
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:501
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1243
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:80
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
label nFaces() const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
label nCells() const
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
static autoPtr< renumberMethod > New(const dictionary &renumberDict)
Return a reference to the selected renumbering method.
label nTotalCells() const
Return total number of cells in decomposed mesh.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
const cellList & cells() const
Field reading functions for post-processing utilities.
const dimensionSet dimless
fvMesh & mesh
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual const labelList & faceNeighbour() const =0
Face face-neighbour addressing.
Abstract base class for renumbering.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:882
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
void Swap(T &a, T &b)
Definition: Swap.H:43
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Cuthill-McKee renumbering.
void sort(UList< T > &)
Definition: UList.C:115
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
static autoPtr< decompositionMethod > NewDecomposer(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1435
static const label labelMax
Definition: label.H:62
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Foam::polyBoundaryMesh.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
ListType reorder(const labelUList &oldToNew, const ListType &)
Reorder the elements (indices, not values) of a list.
static const char nl
Definition: Ostream.H:260
void resetPrimitives(pointField &&points, faceList &&faces, labelList &&owner, labelList &&neighbour, const labelList &patchSizes, const labelList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:672
objects
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Internal & ref()
Return a reference to the dimensioned internal field.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A collection of cell labels.
Definition: cellSet.H:48
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:495
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
A List with indirect addressing.
Definition: fvMatrix.H:106
Direct mesh changes based on v1.3 polyTopoChange syntax.
messageStream Info
virtual const labelList & faceOwner() const =0
Face face-owner addressing.
label nPoints() const
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
A class for managing temporary objects.
Definition: PtrList.H:53
bool rm(const fileName &)
Remove a file, returning true if successful otherwise false.
Definition: POSIX.C:1017
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864