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-2024 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"
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 using namespace Foam;
55 
56 void writeCellLabels
57 (
58  const fvMesh& mesh,
59  const word& name,
60  const labelList& elems
61 )
62 {
64  (
65  IOobject
66  (
67  name,
68  mesh.time().name(),
69  mesh,
72  ),
73  mesh,
74  dimless,
75  scalarField(scalarList(elems))
76  );
77 
78  fld.write();
79 }
80 
81 
82 // Calculate band of matrix
83 label getBand(const labelList& owner, const labelList& neighbour)
84 {
85  label band = 0;
86 
87  forAll(neighbour, facei)
88  {
89  label diff = neighbour[facei] - owner[facei];
90 
91  if (diff > band)
92  {
93  band = diff;
94  }
95  }
96  return band;
97 }
98 
99 
100 // Calculate band of matrix
101 void getBand
102 (
103  const bool calculateIntersect,
104  const label nCells,
105  const labelList& owner,
106  const labelList& neighbour,
107  label& bandwidth,
108  scalar& profile, // scalar to avoid overflow
109  scalar& sumSqrIntersect // scalar to avoid overflow
110 )
111 {
112  labelList cellBandwidth(nCells, 0);
113  scalarField nIntersect(nCells, 0.0);
114 
115  forAll(neighbour, facei)
116  {
117  label own = owner[facei];
118  label nei = neighbour[facei];
119 
120  // Note: mag not necessary for correct (upper-triangular) ordering.
121  label diff = nei-own;
122  cellBandwidth[nei] = max(cellBandwidth[nei], diff);
123  }
124 
125  bandwidth = max(cellBandwidth);
126 
127  // Do not use field algebra because of conversion label to scalar
128  profile = 0.0;
129  forAll(cellBandwidth, celli)
130  {
131  profile += 1.0*cellBandwidth[celli];
132  }
133 
134  sumSqrIntersect = 0.0;
135  if (calculateIntersect)
136  {
137  forAll(nIntersect, celli)
138  {
139  for (label colI = celli-cellBandwidth[celli]; colI <= celli; colI++)
140  {
141  nIntersect[colI] += 1.0;
142  }
143  }
144 
145  sumSqrIntersect = sum(Foam::sqr(nIntersect));
146  }
147 }
148 
149 
150 // Determine upper-triangular face order
151 labelList getFaceOrder
152 (
153  const primitiveMesh& mesh,
154  const labelList& cellOrder // New to old cell
155 )
156 {
157  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
158 
159  labelList oldToNewFace(mesh.nFaces(), -1);
160 
161  label newFacei = 0;
162 
163  labelList nbr;
165 
166  forAll(cellOrder, newCelli)
167  {
168  label oldCelli = cellOrder[newCelli];
169 
170  const cell& cFaces = mesh.cells()[oldCelli];
171 
172  // Neighbouring cells
173  nbr.setSize(cFaces.size());
174 
175  forAll(cFaces, i)
176  {
177  label facei = cFaces[i];
178 
179  if (mesh.isInternalFace(facei))
180  {
181  // Internal face. Get cell on other side.
182  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
183  if (nbrCelli == newCelli)
184  {
185  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
186  }
187 
188  if (newCelli < nbrCelli)
189  {
190  // Celli is master
191  nbr[i] = nbrCelli;
192  }
193  else
194  {
195  // nbrCell is master. Let it handle this face.
196  nbr[i] = -1;
197  }
198  }
199  else
200  {
201  // External face. Do later.
202  nbr[i] = -1;
203  }
204  }
205 
206  order.setSize(nbr.size());
207  sortedOrder(nbr, order);
208 
209  forAll(order, i)
210  {
211  label index = order[i];
212  if (nbr[index] != -1)
213  {
214  oldToNewFace[cFaces[index]] = newFacei++;
215  }
216  }
217  }
218 
219  // Leave patch faces intact.
220  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
221  {
222  oldToNewFace[facei] = facei;
223  }
224 
225 
226  // Check done all faces.
227  forAll(oldToNewFace, facei)
228  {
229  if (oldToNewFace[facei] == -1)
230  {
232  << "Did not determine new position" << " for face " << facei
233  << abort(FatalError);
234  }
235  }
236 
237  return invert(mesh.nFaces(), oldToNewFace);
238 }
239 
240 
241 // Determine face order such that inside region faces are sorted
242 // upper-triangular but in between region faces are handled like boundary faces.
243 labelList getRegionFaceOrder
244 (
245  const primitiveMesh& mesh,
246  const labelList& cellOrder, // New to old cell
247  const labelList& cellToRegion // Old cell to region
248 )
249 {
250  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
251 
252  labelList oldToNewFace(mesh.nFaces(), -1);
253 
254  label newFacei = 0;
255 
256  label prevRegion = -1;
257 
258  forAll(cellOrder, newCelli)
259  {
260  label oldCelli = cellOrder[newCelli];
261 
262  if (cellToRegion[oldCelli] != prevRegion)
263  {
264  prevRegion = cellToRegion[oldCelli];
265  }
266 
267  const cell& cFaces = mesh.cells()[oldCelli];
268 
269  SortableList<label> nbr(cFaces.size());
270 
271  forAll(cFaces, i)
272  {
273  label facei = cFaces[i];
274 
275  if (mesh.isInternalFace(facei))
276  {
277  // Internal face. Get cell on other side.
278  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
279  if (nbrCelli == newCelli)
280  {
281  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
282  }
283 
284  if (cellToRegion[oldCelli] != cellToRegion[cellOrder[nbrCelli]])
285  {
286  // Treat like external face. Do later.
287  nbr[i] = -1;
288  }
289  else if (newCelli < nbrCelli)
290  {
291  // Celli is master
292  nbr[i] = nbrCelli;
293  }
294  else
295  {
296  // nbrCell is master. Let it handle this face.
297  nbr[i] = -1;
298  }
299  }
300  else
301  {
302  // External face. Do later.
303  nbr[i] = -1;
304  }
305  }
306 
307  nbr.sort();
308 
309  forAll(nbr, i)
310  {
311  if (nbr[i] != -1)
312  {
313  oldToNewFace[cFaces[nbr.indices()[i]]] = newFacei++;
314  }
315  }
316  }
317 
318  // Do region interfaces
319  label nRegions = max(cellToRegion)+1;
320  {
321  // Sort in increasing region
322  SortableList<label> sortKey(mesh.nFaces(), labelMax);
323 
324  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
325  {
326  label ownRegion = cellToRegion[mesh.faceOwner()[facei]];
327  label neiRegion = cellToRegion[mesh.faceNeighbour()[facei]];
328 
329  if (ownRegion != neiRegion)
330  {
331  sortKey[facei] =
332  min(ownRegion, neiRegion)*nRegions
333  +max(ownRegion, neiRegion);
334  }
335  }
336  sortKey.sort();
337 
338  // Extract.
339  label prevKey = -1;
340  forAll(sortKey, i)
341  {
342  label key = sortKey[i];
343 
344  if (key == labelMax)
345  {
346  break;
347  }
348 
349  if (prevKey != key)
350  {
351  prevKey = key;
352  }
353 
354  oldToNewFace[sortKey.indices()[i]] = newFacei++;
355  }
356  }
357 
358  // Leave patch faces intact.
359  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
360  {
361  oldToNewFace[facei] = facei;
362  }
363 
364 
365  // Check done all faces.
366  forAll(oldToNewFace, facei)
367  {
368  if (oldToNewFace[facei] == -1)
369  {
371  << "Did not determine new position"
372  << " for face " << facei
373  << abort(FatalError);
374  }
375  }
376 
377  return invert(mesh.nFaces(), oldToNewFace);
378 }
379 
380 
381 // cellOrder: old cell for every new cell
382 // faceOrder: old face for every new face. Ordering of boundary faces not
383 // changed.
384 autoPtr<polyTopoChangeMap> reorderMesh
385 (
386  polyMesh& mesh,
387  labelList& cellOrder,
388  labelList& faceOrder
389 )
390 {
391  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
392  labelList reverseFaceOrder(invert(faceOrder.size(), faceOrder));
393 
394  faceList newFaces(reorder(reverseFaceOrder, mesh.faces()));
395  labelList newOwner
396  (
397  renumber
398  (
399  reverseCellOrder,
400  reorder(reverseFaceOrder, mesh.faceOwner())
401  )
402  );
403  labelList newNeighbour
404  (
405  renumber
406  (
407  reverseCellOrder,
408  reorder(reverseFaceOrder, mesh.faceNeighbour())
409  )
410  );
411 
412  // Check if any faces need swapping.
413  labelHashSet flipFaceFlux(newOwner.size());
414  forAll(newNeighbour, facei)
415  {
416  label own = newOwner[facei];
417  label nei = newNeighbour[facei];
418 
419  if (nei < own)
420  {
421  newFaces[facei].flip();
422  Swap(newOwner[facei], newNeighbour[facei]);
423  flipFaceFlux.insert(facei);
424  }
425  }
426 
427  const polyBoundaryMesh& patches = mesh.boundaryMesh();
428  labelList patchSizes(patches.size());
429  labelList patchStarts(patches.size());
430  labelList oldPatchNMeshPoints(patches.size());
431  labelListList patchPointMap(patches.size());
432 
434  {
435  patchSizes[patchi] = patches[patchi].size();
436  patchStarts[patchi] = patches[patchi].start();
437  oldPatchNMeshPoints[patchi] = patches[patchi].nPoints();
438  patchPointMap[patchi] = identityMap(patches[patchi].nPoints());
439  }
440 
441  mesh.resetPrimitives
442  (
443  NullObjectMove<pointField>(),
444  move(newFaces),
445  move(newOwner),
446  move(newNeighbour),
447  patchSizes,
448  patchStarts,
449  true
450  );
451 
453  (
455  (
456  mesh, // const polyMesh& mesh,
457  mesh.nPoints(), // nOldPoints,
458  mesh.nFaces(), // nOldFaces,
459  mesh.nCells(), // nOldCells,
460  identityMap(mesh.nPoints()), // pointMap,
461  List<objectMap>(0), // pointsFromPoints,
462  move(faceOrder), // faceMap,
463  List<objectMap>(0), // facesFromFaces,
464  move(cellOrder), // cellMap,
465  List<objectMap>(0), // cellsFromCells,
466  identityMap(mesh.nPoints()), // reversePointMap,
467  move(reverseFaceOrder), // reverseFaceMap,
468  move(reverseCellOrder), // reverseCellMap,
469  move(flipFaceFlux), // flipFaceFlux,
470  move(patchPointMap), // patchPointMap,
471  move(patchSizes), // oldPatchSizes
472  move(patchStarts), // oldPatchStarts,
473  move(oldPatchNMeshPoints), // oldPatchNMeshPoints
474  autoPtr<scalarField>() // oldCellVolumes
475  )
476  );
477 }
478 
479 
480 // Return new to old cell numbering
481 labelList regionRenumber
482 (
483  const renumberMethod& method,
484  const fvMesh& mesh,
485  const labelList& cellToRegion
486 )
487 {
488  Info<< "Determining cell order:" << endl;
489 
490  labelList cellOrder(cellToRegion.size());
491 
492  label nRegions = max(cellToRegion)+1;
493 
494  labelListList regionToCells(invertOneToMany(nRegions, cellToRegion));
495 
496  label celli = 0;
497 
498  forAll(regionToCells, regionI)
499  {
500  Info<< " region " << regionI << " starts at " << celli << endl;
501 
502  // Make sure no parallel comms
503  bool oldParRun = UPstream::parRun();
504  UPstream::parRun() = false;
505 
506  // Per region do a reordering.
507  fvMeshSubset subsetter(mesh);
508  subsetter.setLargeCellSubset(cellToRegion, regionI);
509 
510  const fvMesh& subMesh = subsetter.subMesh();
511 
512  labelList subCellOrder = method.renumber
513  (
514  subMesh,
515  subMesh.cellCentres()
516  );
517 
518  // Restore state
519  UPstream::parRun() = oldParRun;
520 
521  const labelList& cellMap = subsetter.cellMap();
522 
523  forAll(subCellOrder, i)
524  {
525  cellOrder[celli++] = cellMap[subCellOrder[i]];
526  }
527  }
528  Info<< endl;
529 
530  return cellOrder;
531 }
532 
533 
534 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
535 
536 int main(int argc, char *argv[])
537 {
539  (
540  "Renumber mesh to minimise bandwidth"
541  );
542 
543  #include "addRegionOption.H"
544  #include "addOverwriteOption.H"
546  #include "addDictOption.H"
548  (
549  "frontWidth",
550  "calculate the rms of the frontwidth"
551  );
553  (
554  "noFields",
555  "do not update fields"
556  );
557 
558  #include "setRootCase.H"
560 
561  // Get times list
562  const instantList Times = timeSelector::select0(runTime, args);
563 
565  const word oldInstance = mesh.pointsInstance();
566 
567  const bool readDict = args.optionFound("dict");
568  const bool doFrontWidth = args.optionFound("frontWidth");
569  const bool overwrite = args.optionFound("overwrite");
570  const bool fields = !args.optionFound("noFields");
571 
572  label band;
573  scalar profile;
574  scalar sumSqrIntersect;
575  getBand
576  (
577  doFrontWidth,
578  mesh.nCells(),
579  mesh.faceOwner(),
580  mesh.faceNeighbour(),
581  band,
582  profile,
583  sumSqrIntersect
584  );
585 
586  reduce(band, maxOp<label>());
587  reduce(profile, sumOp<scalar>());
588  scalar rmsFrontwidth = Foam::sqrt
589  (
591  (
592  sumSqrIntersect,
593  sumOp<scalar>()
594  )/mesh.globalData().nTotalCells()
595  );
596 
597  Info<< "Mesh size: " << mesh.globalData().nTotalCells() << nl
598  << "Before renumbering :" << nl
599  << " band : " << band << nl
600  << " profile : " << profile << nl;
601 
602  if (doFrontWidth)
603  {
604  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
605  }
606 
607  Info<< endl;
608 
609  bool sortCoupledFaceCells = false;
610  bool writeMaps = false;
611  bool orderPoints = false;
612  label blockSize = 0;
613  bool renumberSets = true;
614 
615  // Construct renumberMethod
616  autoPtr<dictionary> renumberDictPtr;
617  autoPtr<renumberMethod> renumberPtr;
618 
619  if (readDict)
620  {
621  renumberDictPtr.reset
622  (
623  new dictionary(systemDict("renumberMeshDict", args, mesh))
624  );
625  const dictionary& renumberDict = renumberDictPtr();
626 
627  renumberPtr = renumberMethod::New(renumberDict);
628 
629  sortCoupledFaceCells = renumberDict.lookupOrDefault
630  (
631  "sortCoupledFaceCells",
632  false
633  );
634  if (sortCoupledFaceCells)
635  {
636  Info<< "Sorting cells on coupled boundaries to be last." << nl
637  << endl;
638  }
639 
640  blockSize = renumberDict.lookupOrDefault("blockSize", 0);
641  if (blockSize > 0)
642  {
643  Info<< "Ordering cells into regions of size " << blockSize
644  << " (using decomposition);"
645  << " ordering faces into region-internal and region-external."
646  << nl << endl;
647 
648  if (blockSize < 0 || blockSize >= mesh.nCells())
649  {
651  << "Block size " << blockSize
652  << " should be positive integer"
653  << " and less than the number of cells in the mesh."
654  << exit(FatalError);
655  }
656  }
657 
658  orderPoints = renumberDict.lookupOrDefault("orderPoints", false);
659  if (orderPoints)
660  {
661  Info<< "Ordering points into internal and boundary points." << nl
662  << endl;
663  }
664 
665  renumberDict.lookup("writeMaps") >> writeMaps;
666  if (writeMaps)
667  {
668  Info<< "Writing renumber maps (new to old) to polyMesh." << nl
669  << endl;
670  }
671 
672  renumberSets = renumberDict.lookupOrDefault("renumberSets", true);
673  }
674  else
675  {
676  Info<< "Using default renumberMethod." << nl << endl;
677  dictionary renumberDict;
678  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
679  }
680 
681  Info<< "Selecting renumberMethod " << renumberPtr().type() << nl << endl;
682 
683 
684 
685  // Read parallel reconstruct maps
686  labelIOList cellProcAddressing
687  (
688  IOobject
689  (
690  "cellProcAddressing",
691  mesh.facesInstance(),
693  mesh,
695  ),
696  labelList(0)
697  );
698 
699  labelIOList faceProcAddressing
700  (
701  IOobject
702  (
703  "faceProcAddressing",
704  mesh.facesInstance(),
706  mesh,
708  ),
709  labelList(0)
710  );
711  labelIOList pointProcAddressing
712  (
713  IOobject
714  (
715  "pointProcAddressing",
716  mesh.pointsInstance(),
718  mesh,
720  ),
721  labelList(0)
722  );
723 
724 
725  // Read objects in time directory
726  IOobjectList objects(mesh, runTime.name());
727 
728  if (fields) Info<< "Reading geometric fields" << nl << endl;
729 
730  #include "readVolFields.H"
731  #include "readSurfaceFields.H"
732  #include "readPointFields.H"
733 
734  Info<< endl;
735 
736  // From renumbering:
737  // - from new cell/face back to original cell/face
738  labelList cellOrder;
739  labelList faceOrder;
740 
741  if (blockSize > 0)
742  {
743  // Renumbering in two phases. Should be done in one so mapping of
744  // fields is done correctly!
745 
746  label nBlocks = mesh.nCells()/blockSize;
747  Info<< "nBlocks = " << nBlocks << endl;
748 
749  // Read decompositionMethod dictionary
750  dictionary decomposeDict(renumberDictPtr().subDict("blockCoeffs"));
751  decomposeDict.set("numberOfSubdomains", nBlocks);
752 
753  bool oldParRun = UPstream::parRun();
754  UPstream::parRun() = false;
755 
756  autoPtr<decompositionMethod> decomposePtr
757  (
759  );
760 
761  labelList cellToRegion
762  (
763  decomposePtr().decompose
764  (
765  mesh,
766  mesh.cellCentres()
767  )
768  );
769 
770  // Restore state
771  UPstream::parRun() = oldParRun;
772 
773  // For debugging: write out region
774  writeCellLabels
775  (
776  mesh,
777  "cellProc",
778  cellToRegion
779  );
780  Info<< nl << "Written decomposition as volScalarField::Internal to "
781  << "cellProc for use in postprocessing."
782  << nl << endl;
783 
784  cellOrder = regionRenumber(renumberPtr(), mesh, cellToRegion);
785 
786  // Determine new to old face order with new cell numbering
787  faceOrder = getRegionFaceOrder
788  (
789  mesh,
790  cellOrder,
791  cellToRegion
792  );
793  }
794  else
795  {
796  // Determines sorted back to original cell ordering
797  cellOrder = renumberPtr().renumber
798  (
799  mesh,
800  mesh.cellCentres()
801  );
802 
803  if (sortCoupledFaceCells)
804  {
805  // Change order so all coupled patch faceCells are at the end.
806  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
807 
808  // Collect all boundary cells on coupled patches
809  label nBndCells = 0;
810  forAll(pbm, patchi)
811  {
812  if (pbm[patchi].coupled())
813  {
814  nBndCells += pbm[patchi].size();
815  }
816  }
817 
818  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
819 
820  labelList bndCells(nBndCells);
821  labelList bndCellMap(nBndCells);
822  nBndCells = 0;
823  forAll(pbm, patchi)
824  {
825  if (pbm[patchi].coupled())
826  {
827  const labelUList& faceCells = pbm[patchi].faceCells();
828  forAll(faceCells, i)
829  {
830  label celli = faceCells[i];
831 
832  if (reverseCellOrder[celli] != -1)
833  {
834  bndCells[nBndCells] = celli;
835  bndCellMap[nBndCells++] = reverseCellOrder[celli];
836  reverseCellOrder[celli] = -1;
837  }
838  }
839  }
840  }
841  bndCells.setSize(nBndCells);
842  bndCellMap.setSize(nBndCells);
843 
844  // Sort
846  sortedOrder(bndCellMap, order);
847 
848  // Redo newReverseCellOrder
849  labelList newReverseCellOrder(mesh.nCells(), -1);
850 
851  label sortedI = mesh.nCells();
852  forAllReverse(order, i)
853  {
854  label origCelli = bndCells[order[i]];
855  newReverseCellOrder[origCelli] = --sortedI;
856  }
857 
858  Info<< "Ordered all " << nBndCells << " cells with a coupled face"
859  << " to the end of the cell list, starting at " << sortedI
860  << endl;
861 
862  // Compact
863  sortedI = 0;
864  forAll(cellOrder, newCelli)
865  {
866  label origCelli = cellOrder[newCelli];
867  if (newReverseCellOrder[origCelli] == -1)
868  {
869  newReverseCellOrder[origCelli] = sortedI++;
870  }
871  }
872 
873  // Update sorted back to original (unsorted) map
874  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
875  }
876 
877 
878  // Determine new to old face order with new cell numbering
879  faceOrder = getFaceOrder
880  (
881  mesh,
882  cellOrder // New to old cell
883  );
884  }
885 
886 
887  if (!overwrite)
888  {
889  runTime++;
890  }
891 
892 
893  // Change the mesh.
894  autoPtr<polyTopoChangeMap> map = reorderMesh(mesh, cellOrder, faceOrder);
895 
896 
897  if (orderPoints)
898  {
899  polyTopoChange meshMod(mesh);
900  autoPtr<polyTopoChangeMap> pointOrderMap = meshMod.changeMesh
901  (
902  mesh,
903  true, // syncParallel
904  false, // orderCells
905  orderPoints // orderPoints
906  );
907 
908  // Combine point reordering into map.
909  const_cast<labelList&>(map().pointMap()) = UIndirectList<label>
910  (
911  map().pointMap(),
912  pointOrderMap().pointMap()
913  )();
914 
916  (
917  pointOrderMap().reversePointMap(),
918  const_cast<labelList&>(map().reversePointMap())
919  );
920  }
921 
922 
923  // Update fields
924  mesh.topoChange(map);
925 
926  // Update proc maps
927  if
928  (
929  cellProcAddressing.headerOk()
930  && cellProcAddressing.size() == mesh.nCells()
931  )
932  {
933  Info<< "Renumbering processor cell decomposition map "
934  << cellProcAddressing.name() << endl;
935 
936  cellProcAddressing = labelList
937  (
938  UIndirectList<label>(cellProcAddressing, map().cellMap())
939  );
940  }
941  if
942  (
943  faceProcAddressing.headerOk()
944  && faceProcAddressing.size() == mesh.nFaces()
945  )
946  {
947  Info<< "Renumbering processor face decomposition map "
948  << faceProcAddressing.name() << endl;
949 
950  faceProcAddressing = labelList
951  (
952  UIndirectList<label>(faceProcAddressing, map().faceMap())
953  );
954 
955  // Detect any flips.
956  const labelHashSet& fff = map().flipFaceFlux();
957  forAllConstIter(labelHashSet, fff, iter)
958  {
959  label facei = iter.key();
960  label masterFacei = faceProcAddressing[facei];
961 
962  faceProcAddressing[facei] = -masterFacei;
963 
964  if (masterFacei == 0)
965  {
967  << " masterFacei:" << masterFacei << exit(FatalError);
968  }
969  }
970  }
971  if
972  (
973  pointProcAddressing.headerOk()
974  && pointProcAddressing.size() == mesh.nPoints()
975  )
976  {
977  Info<< "Renumbering processor point decomposition map "
978  << pointProcAddressing.name() << endl;
979 
980  pointProcAddressing = labelList
981  (
982  UIndirectList<label>(pointProcAddressing, map().pointMap())
983  );
984  }
985 
986 
987  {
988  label band;
989  scalar profile;
990  scalar sumSqrIntersect;
991  getBand
992  (
993  doFrontWidth,
994  mesh.nCells(),
995  mesh.faceOwner(),
996  mesh.faceNeighbour(),
997  band,
998  profile,
999  sumSqrIntersect
1000  );
1001  reduce(band, maxOp<label>());
1002  reduce(profile, sumOp<scalar>());
1003  scalar rmsFrontwidth = Foam::sqrt
1004  (
1005  returnReduce
1006  (
1007  sumSqrIntersect,
1008  sumOp<scalar>()
1009  )/mesh.globalData().nTotalCells()
1010  );
1011 
1012  Info<< "After renumbering :" << nl
1013  << " band : " << band << nl
1014  << " profile : " << profile << nl;
1015 
1016  if (doFrontWidth)
1017  {
1018 
1019  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1020  }
1021 
1022  Info<< endl;
1023  }
1024 
1025  if (orderPoints)
1026  {
1027  // Force edge calculation (since only reason that points would need to
1028  // be sorted)
1029  (void)mesh.edges();
1030 
1031  label nTotPoints = returnReduce
1032  (
1033  mesh.nPoints(),
1034  sumOp<label>()
1035  );
1036  label nTotIntPoints = returnReduce
1037  (
1038  mesh.nInternalPoints(),
1039  sumOp<label>()
1040  );
1041 
1042  label nTotEdges = returnReduce
1043  (
1044  mesh.nEdges(),
1045  sumOp<label>()
1046  );
1047  label nTotIntEdges = returnReduce
1048  (
1049  mesh.nInternalEdges(),
1050  sumOp<label>()
1051  );
1052  label nTotInt0Edges = returnReduce
1053  (
1054  mesh.nInternal0Edges(),
1055  sumOp<label>()
1056  );
1057  label nTotInt1Edges = returnReduce
1058  (
1059  mesh.nInternal1Edges(),
1060  sumOp<label>()
1061  );
1062 
1063  Info<< "Points:" << nl
1064  << " total : " << nTotPoints << nl
1065  << " internal: " << nTotIntPoints << nl
1066  << " boundary: " << nTotPoints-nTotIntPoints << nl
1067  << "Edges:" << nl
1068  << " total : " << nTotEdges << nl
1069  << " internal: " << nTotIntEdges << nl
1070  << " internal using 0 boundary points: "
1071  << nTotInt0Edges << nl
1072  << " internal using 1 boundary points: "
1073  << nTotInt1Edges-nTotInt0Edges << nl
1074  << " internal using 2 boundary points: "
1075  << nTotIntEdges-nTotInt1Edges << nl
1076  << " boundary: " << nTotEdges-nTotIntEdges << nl
1077  << endl;
1078  }
1079 
1080 
1081  if (overwrite)
1082  {
1083  mesh.setInstance(oldInstance);
1084  }
1085 
1086  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1087 
1088  mesh.write();
1089 
1090  if (cellProcAddressing.headerOk())
1091  {
1092  cellProcAddressing.instance() = mesh.facesInstance();
1093  if (cellProcAddressing.size() == mesh.nCells())
1094  {
1095  cellProcAddressing.write();
1096  }
1097  else
1098  {
1099  // procAddressing file no longer valid. Delete it.
1100  const fileName fName(cellProcAddressing.filePath());
1101  if (fName.size())
1102  {
1103  Info<< "Deleting inconsistent processor cell decomposition"
1104  << " map " << fName << endl;
1105  rm(fName);
1106  }
1107  }
1108  }
1109 
1110  if (faceProcAddressing.headerOk())
1111  {
1112  faceProcAddressing.instance() = mesh.facesInstance();
1113  if (faceProcAddressing.size() == mesh.nFaces())
1114  {
1115  faceProcAddressing.write();
1116  }
1117  else
1118  {
1119  const fileName fName(faceProcAddressing.filePath());
1120  if (fName.size())
1121  {
1122  Info<< "Deleting inconsistent processor face decomposition"
1123  << " map " << fName << endl;
1124  rm(fName);
1125  }
1126  }
1127  }
1128 
1129  if (pointProcAddressing.headerOk())
1130  {
1131  pointProcAddressing.instance() = mesh.facesInstance();
1132  if (pointProcAddressing.size() == mesh.nPoints())
1133  {
1134  pointProcAddressing.write();
1135  }
1136  else
1137  {
1138  const fileName fName(pointProcAddressing.filePath());
1139  if (fName.size())
1140  {
1141  Info<< "Deleting inconsistent processor point decomposition"
1142  << " map " << fName << endl;
1143  rm(fName);
1144  }
1145  }
1146  }
1147 
1148  if (writeMaps)
1149  {
1150  // For debugging: write out region
1151  writeCellLabels
1152  (
1153  mesh,
1154  "origCellID",
1155  map().cellMap()
1156  );
1157 
1158  writeCellLabels
1159  (
1160  mesh,
1161  "cellID",
1162  identityMap(mesh.nCells())
1163  );
1164 
1165  Info<< nl << "Written current cellID and origCellID as volScalarField"
1166  << " for use in postprocessing."
1167  << nl << endl;
1168 
1169  labelIOList
1170  (
1171  IOobject
1172  (
1173  "cellMap",
1174  mesh.facesInstance(),
1176  mesh,
1179  false
1180  ),
1181  map().cellMap()
1182  ).write();
1183 
1184  labelIOList
1185  (
1186  IOobject
1187  (
1188  "faceMap",
1189  mesh.facesInstance(),
1191  mesh,
1194  false
1195  ),
1196  map().faceMap()
1197  ).write();
1198 
1199  labelIOList
1200  (
1201  IOobject
1202  (
1203  "pointMap",
1204  mesh.facesInstance(),
1206  mesh,
1209  false
1210  ),
1211  map().pointMap()
1212  ).write();
1213  }
1214 
1215 
1216  // Renumber sets if required
1217  if (renumberSets)
1218  {
1219  Info<< endl;
1220 
1221  // Read sets
1222  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1223 
1224  {
1225  IOobjectList cSets(objects.lookupClass(cellSet::typeName));
1226  if (cSets.size())
1227  {
1228  Info<< "Renumbering cellSets:" << endl;
1229  forAllConstIter(IOobjectList, cSets, iter)
1230  {
1231  cellSet cs(*iter());
1232  Info<< " " << cs.name() << endl;
1233  cs.topoChange(map());
1234  cs.instance() = mesh.facesInstance();
1235  cs.write();
1236  }
1237  }
1238  }
1239 
1240  {
1241  IOobjectList fSets(objects.lookupClass(faceSet::typeName));
1242  if (fSets.size())
1243  {
1244  Info<< "Renumbering faceSets:" << endl;
1245  forAllConstIter(IOobjectList, fSets, iter)
1246  {
1247  faceSet fs(*iter());
1248  Info<< " " << fs.name() << endl;
1249  fs.topoChange(map());
1250  fs.instance() = mesh.facesInstance();
1251  fs.write();
1252  }
1253  }
1254  }
1255 
1256  {
1257  IOobjectList pSets(objects.lookupClass(pointSet::typeName));
1258  if (pSets.size())
1259  {
1260  Info<< "Renumbering pointSets:" << endl;
1261  forAllConstIter(IOobjectList, pSets, iter)
1262  {
1263  pointSet ps(*iter());
1264  Info<< " " << ps.name() << endl;
1265  ps.topoChange(map());
1266  ps.instance() = mesh.facesInstance();
1267  ps.write();
1268  }
1269  }
1270  }
1271  }
1272 
1273  Info<< "End\n" << endl;
1274 
1275  return 0;
1276 }
1277 
1278 
1279 // ************************************************************************* //
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
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:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const word & name() const
Return const reference to name.
A list of face labels.
Definition: faceSet.H:51
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:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1737
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1264
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:971
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
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:703
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1515
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:965
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1345
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:270
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
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:334
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:228
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:257
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
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)
dimensionedScalar sqrt(const dimensionedScalar &ds)
ListType reorder(const label size, const typename ListType::value_type &defaultValue, const labelUList &oldToNew, const ListType &lst)
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:409
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
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:266
objects
Foam::argList args(argc, argv)
Foam::surfaceFields.