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