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