renumberMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  renumberMesh
29 
30 Group
31  grpMeshManipulationUtilities
32 
33 Description
34  Renumbers the cell list in order to reduce the bandwidth, reading and
35  renumbering all fields from all the time directories.
36 
37  By default uses bandCompression (CuthillMcKee) but will
38  read system/renumberMeshDict if -dict option is present
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #include "argList.H"
43 #include "IOobjectList.H"
44 #include "fvMesh.H"
45 #include "polyTopoChange.H"
46 #include "ReadFields.H"
47 #include "volFields.H"
48 #include "surfaceFields.H"
49 #include "SortableList.H"
50 #include "decompositionMethod.H"
51 #include "renumberMethod.H"
53 #include "CuthillMcKeeRenumber.H"
54 #include "fvMeshSubset.H"
55 #include "cellSet.H"
56 #include "faceSet.H"
57 #include "pointSet.H"
58 #include "processorMeshes.H"
59 #include "hexRef8Data.H"
60 #include "regionProperties.H"
61 #include "polyMeshTools.H"
62 
63 #ifdef HAVE_ZOLTAN
64  #include "zoltanRenumber.H"
65 #endif
66 
67 
68 using namespace Foam;
69 
70 
71 // Create named field from labelList for post-processing
72 tmp<volScalarField> createScalarField
73 (
74  const fvMesh& mesh,
75  const word& name,
76  const labelList& elems
77 )
78 {
80  (
81  new volScalarField
82  (
83  IOobject
84  (
85  name,
86  mesh.time().timeName(),
87  mesh,
90  false
91  ),
92  mesh,
94  zeroGradientFvPatchScalarField::typeName
95  )
96  );
97  volScalarField& fld = tfld.ref();
98 
99  forAll(fld, celli)
100  {
101  fld[celli] = elems[celli];
102  }
103 
104  return tfld;
105 }
106 
107 
108 // Calculate band of matrix
109 label getBand(const labelList& owner, const labelList& neighbour)
110 {
111  label band = 0;
112 
113  forAll(neighbour, facei)
114  {
115  label diff = neighbour[facei] - owner[facei];
116 
117  if (diff > band)
118  {
119  band = diff;
120  }
121  }
122  return band;
123 }
124 
125 
126 // Calculate band of matrix
127 void getBand
128 (
129  const bool calculateIntersect,
130  const label nCells,
131  const labelList& owner,
132  const labelList& neighbour,
133  label& bandwidth,
134  scalar& profile, // scalar to avoid overflow
135  scalar& sumSqrIntersect // scalar to avoid overflow
136 )
137 {
138  labelList cellBandwidth(nCells, Zero);
139  scalarField nIntersect(nCells, Zero);
140 
141  forAll(neighbour, facei)
142  {
143  label own = owner[facei];
144  label nei = neighbour[facei];
145 
146  // Note: mag not necessary for correct (upper-triangular) ordering.
147  label diff = nei-own;
148  cellBandwidth[nei] = max(cellBandwidth[nei], diff);
149  }
150 
151  bandwidth = max(cellBandwidth);
152 
153  // Do not use field algebra because of conversion label to scalar
154  profile = 0.0;
155  forAll(cellBandwidth, celli)
156  {
157  profile += 1.0*cellBandwidth[celli];
158  }
159 
160  sumSqrIntersect = 0.0;
161  if (calculateIntersect)
162  {
163  forAll(nIntersect, celli)
164  {
165  for (label colI = celli-cellBandwidth[celli]; colI <= celli; colI++)
166  {
167  nIntersect[colI] += 1.0;
168  }
169  }
170 
171  sumSqrIntersect = sum(Foam::sqr(nIntersect));
172  }
173 }
174 
175 
176 // Determine upper-triangular face order
177 labelList getFaceOrder
178 (
179  const primitiveMesh& mesh,
180  const labelList& cellOrder // New to old cell
181 )
182 {
183  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
184 
185  labelList oldToNewFace(mesh.nFaces(), -1);
186 
187  label newFacei = 0;
188 
189  labelList nbr;
190  labelList order;
191 
192  forAll(cellOrder, newCelli)
193  {
194  label oldCelli = cellOrder[newCelli];
195 
196  const cell& cFaces = mesh.cells()[oldCelli];
197 
198  // Neighbouring cells
199  nbr.setSize(cFaces.size());
200 
201  forAll(cFaces, i)
202  {
203  label facei = cFaces[i];
204 
205  if (mesh.isInternalFace(facei))
206  {
207  // Internal face. Get cell on other side.
208  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
209  if (nbrCelli == newCelli)
210  {
211  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
212  }
213 
214  if (newCelli < nbrCelli)
215  {
216  // Celli is master
217  nbr[i] = nbrCelli;
218  }
219  else
220  {
221  // nbrCell is master. Let it handle this face.
222  nbr[i] = -1;
223  }
224  }
225  else
226  {
227  // External face. Do later.
228  nbr[i] = -1;
229  }
230  }
231 
232  sortedOrder(nbr, order);
233 
234  for (const label index : order)
235  {
236  if (nbr[index] != -1)
237  {
238  oldToNewFace[cFaces[index]] = newFacei++;
239  }
240  }
241  }
242 
243  // Leave patch faces intact.
244  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
245  {
246  oldToNewFace[facei] = facei;
247  }
248 
249 
250  // Check done all faces.
251  forAll(oldToNewFace, facei)
252  {
253  if (oldToNewFace[facei] == -1)
254  {
256  << "Did not determine new position" << " for face " << facei
257  << abort(FatalError);
258  }
259  }
260 
261  return invert(mesh.nFaces(), oldToNewFace);
262 }
263 
264 
265 // Determine face order such that inside region faces are sorted
266 // upper-triangular but inbetween region faces are handled like boundary faces.
267 labelList getRegionFaceOrder
268 (
269  const primitiveMesh& mesh,
270  const labelList& cellOrder, // New to old cell
271  const labelList& cellToRegion // Old cell to region
272 )
273 {
274  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
275 
276  labelList oldToNewFace(mesh.nFaces(), -1);
277 
278  label newFacei = 0;
279 
280  label prevRegion = -1;
281 
282  forAll(cellOrder, newCelli)
283  {
284  label oldCelli = cellOrder[newCelli];
285 
286  if (cellToRegion[oldCelli] != prevRegion)
287  {
288  prevRegion = cellToRegion[oldCelli];
289  }
290 
291  const cell& cFaces = mesh.cells()[oldCelli];
292 
293  SortableList<label> nbr(cFaces.size());
294 
295  forAll(cFaces, i)
296  {
297  label facei = cFaces[i];
298 
299  if (mesh.isInternalFace(facei))
300  {
301  // Internal face. Get cell on other side.
302  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
303  if (nbrCelli == newCelli)
304  {
305  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
306  }
307 
308  if (cellToRegion[oldCelli] != cellToRegion[cellOrder[nbrCelli]])
309  {
310  // Treat like external face. Do later.
311  nbr[i] = -1;
312  }
313  else if (newCelli < nbrCelli)
314  {
315  // Celli is master
316  nbr[i] = nbrCelli;
317  }
318  else
319  {
320  // nbrCell is master. Let it handle this face.
321  nbr[i] = -1;
322  }
323  }
324  else
325  {
326  // External face. Do later.
327  nbr[i] = -1;
328  }
329  }
330 
331  nbr.sort();
332 
333  forAll(nbr, i)
334  {
335  if (nbr[i] != -1)
336  {
337  oldToNewFace[cFaces[nbr.indices()[i]]] = newFacei++;
338  }
339  }
340  }
341 
342  // Do region interfaces
343  label nRegions = max(cellToRegion)+1;
344  {
345  // Sort in increasing region
347 
348  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
349  {
350  label ownRegion = cellToRegion[mesh.faceOwner()[facei]];
351  label neiRegion = cellToRegion[mesh.faceNeighbour()[facei]];
352 
353  if (ownRegion != neiRegion)
354  {
355  sortKey[facei] =
356  min(ownRegion, neiRegion)*nRegions
357  +max(ownRegion, neiRegion);
358  }
359  }
360  sortKey.sort();
361 
362  // Extract.
363  label prevKey = -1;
364  forAll(sortKey, i)
365  {
366  label key = sortKey[i];
367 
368  if (key == labelMax)
369  {
370  break;
371  }
372 
373  if (prevKey != key)
374  {
375  prevKey = key;
376  }
377 
378  oldToNewFace[sortKey.indices()[i]] = newFacei++;
379  }
380  }
381 
382  // Leave patch faces intact.
383  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
384  {
385  oldToNewFace[facei] = facei;
386  }
387 
388 
389  // Check done all faces.
390  forAll(oldToNewFace, facei)
391  {
392  if (oldToNewFace[facei] == -1)
393  {
395  << "Did not determine new position"
396  << " for face " << facei
397  << abort(FatalError);
398  }
399  }
400 
401  return invert(mesh.nFaces(), oldToNewFace);
402 }
403 
404 
405 // cellOrder: old cell for every new cell
406 // faceOrder: old face for every new face. Ordering of boundary faces not
407 // changed.
408 autoPtr<mapPolyMesh> reorderMesh
409 (
410  polyMesh& mesh,
411  const labelList& cellOrder,
412  const labelList& faceOrder
413 )
414 {
415  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
416  labelList reverseFaceOrder(invert(faceOrder.size(), faceOrder));
417 
418  faceList newFaces(reorder(reverseFaceOrder, mesh.faces()));
419  labelList newOwner
420  (
421  renumber
422  (
423  reverseCellOrder,
424  reorder(reverseFaceOrder, mesh.faceOwner())
425  )
426  );
427  labelList newNeighbour
428  (
429  renumber
430  (
431  reverseCellOrder,
432  reorder(reverseFaceOrder, mesh.faceNeighbour())
433  )
434  );
435 
436  // Check if any faces need swapping.
437  labelHashSet flipFaceFlux(newOwner.size());
438  forAll(newNeighbour, facei)
439  {
440  label own = newOwner[facei];
441  label nei = newNeighbour[facei];
442 
443  if (nei < own)
444  {
445  newFaces[facei].flip();
446  std::swap(newOwner[facei], newNeighbour[facei]);
447  flipFaceFlux.insert(facei);
448  }
449  }
450 
452  labelList patchSizes(patches.size());
453  labelList patchStarts(patches.size());
454  labelList oldPatchNMeshPoints(patches.size());
455  labelListList patchPointMap(patches.size());
456 
457  forAll(patches, patchi)
458  {
459  patchSizes[patchi] = patches[patchi].size();
460  patchStarts[patchi] = patches[patchi].start();
461  oldPatchNMeshPoints[patchi] = patches[patchi].nPoints();
462  patchPointMap[patchi] = identity(patches[patchi].nPoints());
463  }
464 
466  (
467  autoPtr<pointField>(), // <- null: leaves points untouched
468  autoPtr<faceList>::New(std::move(newFaces)),
469  autoPtr<labelList>::New(std::move(newOwner)),
470  autoPtr<labelList>::New(std::move(newNeighbour)),
471  patchSizes,
472  patchStarts,
473  true
474  );
475 
476 
477  // Re-do the faceZones
478  {
479  faceZoneMesh& faceZones = mesh.faceZones();
480  faceZones.clearAddressing();
481  forAll(faceZones, zoneI)
482  {
483  faceZone& fZone = faceZones[zoneI];
484  labelList newAddressing(fZone.size());
485  boolList newFlipMap(fZone.size());
486  forAll(fZone, i)
487  {
488  label oldFacei = fZone[i];
489  newAddressing[i] = reverseFaceOrder[oldFacei];
490  if (flipFaceFlux.found(newAddressing[i]))
491  {
492  newFlipMap[i] = !fZone.flipMap()[i];
493  }
494  else
495  {
496  newFlipMap[i] = fZone.flipMap()[i];
497  }
498  }
499  labelList newToOld(sortedOrder(newAddressing));
500  fZone.resetAddressing
501  (
502  labelUIndList(newAddressing, newToOld)(),
503  boolUIndList(newFlipMap, newToOld)()
504  );
505  }
506  }
507  // Re-do the cellZones
508  {
509  cellZoneMesh& cellZones = mesh.cellZones();
510  cellZones.clearAddressing();
511  forAll(cellZones, zoneI)
512  {
513  cellZones[zoneI] = labelUIndList
514  (
515  reverseCellOrder,
516  cellZones[zoneI]
517  )();
518  Foam::sort(cellZones[zoneI]);
519  }
520  }
521 
522 
524  (
525  mesh, // const polyMesh& mesh,
526  mesh.nPoints(), // nOldPoints,
527  mesh.nFaces(), // nOldFaces,
528  mesh.nCells(), // nOldCells,
529  identity(mesh.nPoints()), // pointMap,
530  List<objectMap>(), // pointsFromPoints,
531  faceOrder, // faceMap,
532  List<objectMap>(), // facesFromPoints,
533  List<objectMap>(), // facesFromEdges,
534  List<objectMap>(), // facesFromFaces,
535  cellOrder, // cellMap,
536  List<objectMap>(), // cellsFromPoints,
537  List<objectMap>(), // cellsFromEdges,
538  List<objectMap>(), // cellsFromFaces,
539  List<objectMap>(), // cellsFromCells,
540  identity(mesh.nPoints()), // reversePointMap,
541  reverseFaceOrder, // reverseFaceMap,
542  reverseCellOrder, // reverseCellMap,
543  flipFaceFlux, // flipFaceFlux,
544  patchPointMap, // patchPointMap,
545  labelListList(), // pointZoneMap,
546  labelListList(), // faceZonePointMap,
547  labelListList(), // faceZoneFaceMap,
548  labelListList(), // cellZoneMap,
549  pointField(), // preMotionPoints,
550  patchStarts, // oldPatchStarts,
551  oldPatchNMeshPoints, // oldPatchNMeshPoints
552  autoPtr<scalarField>() // oldCellVolumes
553  );
554 }
555 
556 
557 // Return new to old cell numbering
558 labelList regionRenumber
559 (
560  const renumberMethod& method,
561  const fvMesh& mesh,
562  const labelList& cellToRegion
563 )
564 {
565  Info<< "Determining cell order:" << endl;
566 
567  labelList cellOrder(cellToRegion.size());
568 
569  label nRegions = max(cellToRegion)+1;
570 
571  labelListList regionToCells(invertOneToMany(nRegions, cellToRegion));
572 
573  label celli = 0;
574 
575  forAll(regionToCells, regioni)
576  {
577  Info<< " region " << regioni << " starts at " << celli << endl;
578 
579  // Make sure no parallel comms
580  const bool oldParRun = UPstream::parRun(false);
581 
582  // Per region do a reordering.
583  fvMeshSubset subsetter(mesh, regioni, cellToRegion);
584 
585  const fvMesh& subMesh = subsetter.subMesh();
586 
587  labelList subCellOrder = method.renumber
588  (
589  subMesh,
590  subMesh.cellCentres()
591  );
592 
593  UPstream::parRun(oldParRun); // Restore parallel state
594 
595 
596  const labelList& cellMap = subsetter.cellMap();
597 
598  forAll(subCellOrder, i)
599  {
600  cellOrder[celli++] = cellMap[subCellOrder[i]];
601  }
602  }
603  Info<< endl;
604 
605  return cellOrder;
606 }
607 
608 
609 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
610 
611 int main(int argc, char *argv[])
612 {
614  (
615  "Renumber mesh cells to reduce the bandwidth"
616  );
617 
618  #include "addAllRegionOptions.H"
619  #include "addOverwriteOption.H"
620  #include "addTimeOptions.H"
621 
622  argList::addOption("dict", "file", "Alternative renumberMeshDict");
623 
625  (
626  "frontWidth",
627  "Calculate the rms of the front-width"
628  );
629 
630  argList::noFunctionObjects(); // Never use function objects
631 
632  #include "setRootCase.H"
633  #include "createTime.H"
634  #include "getAllRegionOptions.H"
635 
636  // Force linker to include zoltan symbols. This section is only needed since
637  // Zoltan is a static library
638  #ifdef HAVE_ZOLTAN
639  Info<< "renumberMesh built with zoltan support." << nl << endl;
640  (void)zoltanRenumber::typeName;
641  #endif
642 
643  // Get times list
644  instantList Times = runTime.times();
645 
646  // Set startTime and endTime depending on -time and -latestTime options
647  #include "checkTimeOptions.H"
648 
650 
651  #include "createNamedMeshes.H"
652 
653  const bool readDict = args.found("dict");
654  const bool doFrontWidth = args.found("frontWidth");
655  const bool overwrite = args.found("overwrite");
656 
657 
658  for (fvMesh& mesh : meshes)
659  {
660  // Reset time in case of multiple meshes and not overwrite
662 
663  const word oldInstance = mesh.pointsInstance();
664 
665  label band;
666  scalar profile;
667  scalar sumSqrIntersect;
668  getBand
669  (
670  doFrontWidth,
671  mesh.nCells(),
672  mesh.faceOwner(),
674  band,
675  profile,
676  sumSqrIntersect
677  );
678 
679  reduce(band, maxOp<label>());
680  reduce(profile, sumOp<scalar>());
681  scalar rmsFrontwidth = Foam::sqrt
682  (
684  (
685  sumSqrIntersect,
686  sumOp<scalar>()
688  );
689 
690  Info<< "Mesh " << mesh.name()
691  << " size: " << mesh.globalData().nTotalCells() << nl
692  << "Before renumbering :" << nl
693  << " band : " << band << nl
694  << " profile : " << profile << nl;
695 
696  if (doFrontWidth)
697  {
698  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
699  }
700 
701  Info<< endl;
702 
703  bool sortCoupledFaceCells = false;
704  bool writeMaps = false;
705  bool orderPoints = false;
706  label blockSize = 0;
707 
708  // Construct renumberMethod
709  autoPtr<IOdictionary> renumberDictPtr;
710  autoPtr<renumberMethod> renumberPtr;
711 
712  if (readDict)
713  {
714  const word dictName("renumberMeshDict");
715  #include "setSystemMeshDictionaryIO.H"
716 
717  Info<< "Renumber according to " << dictIO.name() << nl << endl;
718 
719  renumberDictPtr.reset(new IOdictionary(dictIO));
720  const IOdictionary& renumberDict = renumberDictPtr();
721 
722  renumberPtr = renumberMethod::New(renumberDict);
723 
724  sortCoupledFaceCells = renumberDict.getOrDefault
725  (
726  "sortCoupledFaceCells",
727  false
728  );
729  if (sortCoupledFaceCells)
730  {
731  Info<< "Sorting cells on coupled boundaries to be last." << nl
732  << endl;
733  }
734 
735  blockSize = renumberDict.getOrDefault("blockSize", 0);
736  if (blockSize > 0)
737  {
738  Info<< "Ordering cells into regions of size " << blockSize
739  << " (using decomposition);"
740  << " ordering faces into region-internal"
741  << " and region-external."
742  << nl << endl;
743 
744  if (blockSize < 0 || blockSize >= mesh.nCells())
745  {
747  << "Block size " << blockSize
748  << " should be positive integer"
749  << " and less than the number of cells in the mesh."
750  << exit(FatalError);
751  }
752  }
753 
754  orderPoints = renumberDict.getOrDefault("orderPoints", false);
755  if (orderPoints)
756  {
757  Info<< "Ordering points into internal and boundary points."
758  << nl
759  << endl;
760  }
761 
762  renumberDict.readEntry("writeMaps", writeMaps);
763  if (writeMaps)
764  {
765  Info<< "Writing renumber maps (new to old) to polyMesh." << nl
766  << endl;
767  }
768  }
769  else
770  {
771  Info<< "Using default renumberMethod." << nl << endl;
772  dictionary renumberDict;
773  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
774  }
775 
776  Info<< "Selecting renumberMethod " << renumberPtr().type() << nl
777  << endl;
778 
779 
780 
781  // Read parallel reconstruct maps
782  labelIOList cellProcAddressing
783  (
784  IOobject
785  (
786  "cellProcAddressing",
789  mesh,
792  ),
793  labelList(0)
794  );
795 
796  labelIOList faceProcAddressing
797  (
798  IOobject
799  (
800  "faceProcAddressing",
803  mesh,
806  ),
807  labelList(0)
808  );
809  labelIOList pointProcAddressing
810  (
811  IOobject
812  (
813  "pointProcAddressing",
816  mesh,
819  ),
820  labelList(0)
821  );
822  labelIOList boundaryProcAddressing
823  (
824  IOobject
825  (
826  "boundaryProcAddressing",
829  mesh,
832  ),
833  labelList(0)
834  );
835 
836 
837  // Read objects in time directory
838  IOobjectList objects(mesh, runTime.timeName());
839 
840 
841  // Read vol fields.
842 
844  ReadFields(mesh, objects, vsFlds);
845 
847  ReadFields(mesh, objects, vvFlds);
848 
850  ReadFields(mesh, objects, vstFlds);
851 
852  PtrList<volSymmTensorField> vsymtFlds;
853  ReadFields(mesh, objects, vsymtFlds);
854 
856  ReadFields(mesh, objects, vtFlds);
857 
858 
859  // Read surface fields.
860 
862  ReadFields(mesh, objects, ssFlds);
863 
865  ReadFields(mesh, objects, svFlds);
866 
868  ReadFields(mesh, objects, sstFlds);
869 
871  ReadFields(mesh, objects, ssymtFlds);
872 
874  ReadFields(mesh, objects, stFlds);
875 
876 
877  // Read point fields.
878 
880  ReadFields(pointMesh::New(mesh), objects, psFlds);
881 
883  ReadFields(pointMesh::New(mesh), objects, pvFlds);
884 
886  ReadFields(pointMesh::New(mesh), objects, pstFlds);
887 
889  ReadFields(pointMesh::New(mesh), objects, psymtFlds);
890 
892  ReadFields(pointMesh::New(mesh), objects, ptFlds);
893 
894 
895  // Read sets
896  PtrList<cellSet> cellSets;
897  PtrList<faceSet> faceSets;
898  PtrList<pointSet> pointSets;
899  {
900  // Read sets
901  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
902  ReadFields(objects, cellSets);
903  ReadFields(objects, faceSets);
904  ReadFields(objects, pointSets);
905  }
906 
907 
908  Info<< endl;
909 
910  // From renumbering:
911  // - from new cell/face back to original cell/face
912  labelList cellOrder;
913  labelList faceOrder;
914 
915  if (blockSize > 0)
916  {
917  // Renumbering in two phases. Should be done in one so mapping of
918  // fields is done correctly!
919 
920  label nBlocks = mesh.nCells()/blockSize;
921  Info<< "nBlocks = " << nBlocks << endl;
922 
923  // Read decompositionMethod dictionary
924  dictionary decomposeDict(renumberDictPtr().subDict("blockCoeffs"));
925  decomposeDict.set("numberOfSubdomains", nBlocks);
926 
927  const bool oldParRun = UPstream::parRun(false);
928 
930  (
931  decomposeDict
932  );
933 
934  labelList cellToRegion
935  (
936  decomposePtr().decompose
937  (
938  mesh,
939  mesh.cellCentres()
940  )
941  );
942 
943  UPstream::parRun(oldParRun); // Restore parallel state
944 
945 
946  // For debugging: write out region
947  createScalarField
948  (
949  mesh,
950  "cellDist",
951  cellToRegion
952  )().write();
953 
954  Info<< nl << "Written decomposition as volScalarField to "
955  << "cellDist for use in postprocessing."
956  << nl << endl;
957 
958 
959  cellOrder = regionRenumber(renumberPtr(), mesh, cellToRegion);
960 
961  // Determine new to old face order with new cell numbering
962  faceOrder = getRegionFaceOrder
963  (
964  mesh,
965  cellOrder,
966  cellToRegion
967  );
968  }
969  else
970  {
971  // Determines sorted back to original cell ordering
972  cellOrder = renumberPtr().renumber
973  (
974  mesh,
975  mesh.cellCentres()
976  );
977 
978  if (sortCoupledFaceCells)
979  {
980  // Change order so all coupled patch faceCells are at the end.
981  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
982 
983  // Collect all boundary cells on coupled patches
984  label nBndCells = 0;
985  forAll(pbm, patchi)
986  {
987  if (pbm[patchi].coupled())
988  {
989  nBndCells += pbm[patchi].size();
990  }
991  }
992 
993  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
994 
995  labelList bndCells(nBndCells);
996  labelList bndCellMap(nBndCells);
997  nBndCells = 0;
998  forAll(pbm, patchi)
999  {
1000  if (pbm[patchi].coupled())
1001  {
1002  const labelUList& faceCells = pbm[patchi].faceCells();
1003  forAll(faceCells, i)
1004  {
1005  label celli = faceCells[i];
1006 
1007  if (reverseCellOrder[celli] != -1)
1008  {
1009  bndCells[nBndCells] = celli;
1010  bndCellMap[nBndCells++] =
1011  reverseCellOrder[celli];
1012  reverseCellOrder[celli] = -1;
1013  }
1014  }
1015  }
1016  }
1017  bndCells.setSize(nBndCells);
1018  bndCellMap.setSize(nBndCells);
1019 
1020  // Sort
1021  labelList order(sortedOrder(bndCellMap));
1022 
1023  // Redo newReverseCellOrder
1024  labelList newReverseCellOrder(mesh.nCells(), -1);
1025 
1026  label sortedI = mesh.nCells();
1027  forAllReverse(order, i)
1028  {
1029  label origCelli = bndCells[order[i]];
1030  newReverseCellOrder[origCelli] = --sortedI;
1031  }
1032 
1033  Info<< "Ordered all " << nBndCells
1034  << " cells with a coupled face"
1035  << " to the end of the cell list, starting at " << sortedI
1036  << endl;
1037 
1038  // Compact
1039  sortedI = 0;
1040  forAll(cellOrder, newCelli)
1041  {
1042  label origCelli = cellOrder[newCelli];
1043  if (newReverseCellOrder[origCelli] == -1)
1044  {
1045  newReverseCellOrder[origCelli] = sortedI++;
1046  }
1047  }
1048 
1049  // Update sorted back to original (unsorted) map
1050  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
1051  }
1052 
1053 
1054  // Determine new to old face order with new cell numbering
1055  faceOrder = getFaceOrder
1056  (
1057  mesh,
1058  cellOrder // New to old cell
1059  );
1060  }
1061 
1062 
1063  if (!overwrite)
1064  {
1065  ++runTime;
1066  }
1067 
1068 
1069  // Change the mesh.
1070  autoPtr<mapPolyMesh> map = reorderMesh(mesh, cellOrder, faceOrder);
1071 
1072 
1073  if (orderPoints)
1074  {
1075  polyTopoChange meshMod(mesh);
1076  autoPtr<mapPolyMesh> pointOrderMap = meshMod.changeMesh
1077  (
1078  mesh,
1079  false, // inflate
1080  true, // syncParallel
1081  false, // orderCells
1082  orderPoints // orderPoints
1083  );
1084 
1085  // Combine point reordering into map.
1086  const_cast<labelList&>(map().pointMap()) = labelUIndList
1087  (
1088  map().pointMap(),
1089  pointOrderMap().pointMap()
1090  )();
1091 
1093  (
1094  pointOrderMap().reversePointMap(),
1095  const_cast<labelList&>(map().reversePointMap())
1096  );
1097  }
1098 
1099 
1100  // Update fields
1101  mesh.updateMesh(map());
1102 
1103  // Update proc maps
1104  if (cellProcAddressing.headerOk())
1105  {
1106  if (returnReduceAnd(cellProcAddressing.size() == mesh.nCells()))
1107  {
1108  Info<< "Renumbering processor cell decomposition map "
1109  << cellProcAddressing.name() << endl;
1110 
1111  cellProcAddressing = labelList
1112  (
1113  labelUIndList(cellProcAddressing, map().cellMap())
1114  );
1115  }
1116  else
1117  {
1118  Info<< "Not writing inconsistent processor cell decomposition"
1119  << " map " << cellProcAddressing.filePath() << endl;
1120  cellProcAddressing.writeOpt(IOobject::NO_WRITE);
1121  }
1122  }
1123  else
1124  {
1125  cellProcAddressing.writeOpt(IOobject::NO_WRITE);
1126  }
1127 
1128  if (faceProcAddressing.headerOk())
1129  {
1130  if (returnReduceAnd(faceProcAddressing.size() == mesh.nFaces()))
1131  {
1132  Info<< "Renumbering processor face decomposition map "
1133  << faceProcAddressing.name() << endl;
1134 
1135  faceProcAddressing = labelList
1136  (
1137  labelUIndList(faceProcAddressing, map().faceMap())
1138  );
1139 
1140  // Detect any flips.
1141  const labelHashSet& fff = map().flipFaceFlux();
1142  for (const label facei : fff)
1143  {
1144  label masterFacei = faceProcAddressing[facei];
1145 
1146  faceProcAddressing[facei] = -masterFacei;
1147 
1148  if (masterFacei == 0)
1149  {
1151  << " masterFacei:" << masterFacei
1152  << exit(FatalError);
1153  }
1154  }
1155  }
1156  else
1157  {
1158  Info<< "Not writing inconsistent processor face decomposition"
1159  << " map " << faceProcAddressing.filePath() << endl;
1160  faceProcAddressing.writeOpt(IOobject::NO_WRITE);
1161  }
1162  }
1163  else
1164  {
1165  faceProcAddressing.writeOpt(IOobject::NO_WRITE);
1166  }
1167 
1168  if (pointProcAddressing.headerOk())
1169  {
1170  if (returnReduceAnd(pointProcAddressing.size() == mesh.nPoints()))
1171  {
1172  Info<< "Renumbering processor point decomposition map "
1173  << pointProcAddressing.name() << endl;
1174 
1175  pointProcAddressing = labelList
1176  (
1177  labelUIndList(pointProcAddressing, map().pointMap())
1178  );
1179  }
1180  else
1181  {
1182  Info<< "Not writing inconsistent processor point decomposition"
1183  << " map " << pointProcAddressing.filePath() << endl;
1184  pointProcAddressing.writeOpt(IOobject::NO_WRITE);
1185  }
1186  }
1187  else
1188  {
1189  pointProcAddressing.writeOpt(IOobject::NO_WRITE);
1190  }
1191 
1192  if (boundaryProcAddressing.headerOk())
1193  {
1194  if
1195  (
1197  (
1198  boundaryProcAddressing.size() == mesh.boundaryMesh().size()
1199  )
1200  )
1201  {
1202  // No renumbering needed
1203  }
1204  else
1205  {
1206  Info<< "Not writing inconsistent processor patch decomposition"
1207  << " map " << boundaryProcAddressing.filePath() << endl;
1208  boundaryProcAddressing.writeOpt(IOobject::NO_WRITE);
1209  }
1210  }
1211  else
1212  {
1213  boundaryProcAddressing.writeOpt(IOobject::NO_WRITE);
1214  }
1215 
1216 
1217 
1218 
1219  // Move mesh (since morphing might not do this)
1220  if (map().hasMotionPoints())
1221  {
1222  mesh.movePoints(map().preMotionPoints());
1223  }
1224 
1225 
1226  {
1227  label band;
1228  scalar profile;
1229  scalar sumSqrIntersect;
1230  getBand
1231  (
1232  doFrontWidth,
1233  mesh.nCells(),
1234  mesh.faceOwner(),
1235  mesh.faceNeighbour(),
1236  band,
1237  profile,
1238  sumSqrIntersect
1239  );
1240  reduce(band, maxOp<label>());
1241  reduce(profile, sumOp<scalar>());
1242  scalar rmsFrontwidth = Foam::sqrt
1243  (
1244  returnReduce
1245  (
1246  sumSqrIntersect,
1247  sumOp<scalar>()
1249  );
1250 
1251  Info<< "After renumbering :" << nl
1252  << " band : " << band << nl
1253  << " profile : " << profile << nl;
1254 
1255  if (doFrontWidth)
1256  {
1257 
1258  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1259  }
1260 
1261  Info<< endl;
1262  }
1263 
1264  if (orderPoints)
1265  {
1266  // Force edge calculation (since only reason that points would
1267  // need to be sorted)
1268  (void)mesh.edges();
1269 
1270  label nTotPoints = returnReduce
1271  (
1272  mesh.nPoints(),
1273  sumOp<label>()
1274  );
1275  label nTotIntPoints = returnReduce
1276  (
1278  sumOp<label>()
1279  );
1280 
1281  label nTotEdges = returnReduce
1282  (
1283  mesh.nEdges(),
1284  sumOp<label>()
1285  );
1286  label nTotIntEdges = returnReduce
1287  (
1288  mesh.nInternalEdges(),
1289  sumOp<label>()
1290  );
1291  label nTotInt0Edges = returnReduce
1292  (
1294  sumOp<label>()
1295  );
1296  label nTotInt1Edges = returnReduce
1297  (
1299  sumOp<label>()
1300  );
1301 
1302  Info<< "Points:" << nl
1303  << " total : " << nTotPoints << nl
1304  << " internal: " << nTotIntPoints << nl
1305  << " boundary: " << nTotPoints-nTotIntPoints << nl
1306  << "Edges:" << nl
1307  << " total : " << nTotEdges << nl
1308  << " internal: " << nTotIntEdges << nl
1309  << " internal using 0 boundary points: "
1310  << nTotInt0Edges << nl
1311  << " internal using 1 boundary points: "
1312  << nTotInt1Edges-nTotInt0Edges << nl
1313  << " internal using 2 boundary points: "
1314  << nTotIntEdges-nTotInt1Edges << nl
1315  << " boundary: " << nTotEdges-nTotIntEdges << nl
1316  << endl;
1317  }
1318 
1319 
1320  if (overwrite)
1321  {
1322  mesh.setInstance(oldInstance);
1323  }
1324  else
1325  {
1327  }
1328 
1329 
1330  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1331 
1332  // Remove old procAddressing files
1334 
1335  // Update refinement data
1336  hexRef8Data refData
1337  (
1338  IOobject
1339  (
1340  "dummy",
1341  mesh.facesInstance(),
1343  mesh,
1346  false
1347  )
1348  );
1349  refData.updateMesh(map());
1350  refData.write();
1351 
1352  // Update sets
1353  topoSet::updateMesh(mesh.facesInstance(), map(), cellSets);
1354  topoSet::updateMesh(mesh.facesInstance(), map(), faceSets);
1355  topoSet::updateMesh(mesh.facesInstance(), map(), pointSets);
1356 
1357  mesh.write();
1358 
1359  if (writeMaps)
1360  {
1361  // For debugging: write out region
1362  createScalarField
1363  (
1364  mesh,
1365  "origCellID",
1366  map().cellMap()
1367  )().write();
1368 
1369  createScalarField
1370  (
1371  mesh,
1372  "cellID",
1373  identity(mesh.nCells())
1374  )().write();
1375 
1376  Info<< nl
1377  << "Written current cellID and origCellID as volScalarField"
1378  << " for use in postprocessing." << nl << endl;
1379 
1380  labelIOList
1381  (
1382  IOobject
1383  (
1384  "cellMap",
1385  mesh.facesInstance(),
1387  mesh,
1390  false
1391  ),
1392  map().cellMap()
1393  ).write();
1394 
1395  labelIOList
1396  (
1397  IOobject
1398  (
1399  "faceMap",
1400  mesh.facesInstance(),
1402  mesh,
1405  false
1406  ),
1407  map().faceMap()
1408  ).write();
1409 
1410  labelIOList
1411  (
1412  IOobject
1413  (
1414  "pointMap",
1415  mesh.facesInstance(),
1417  mesh,
1420  false
1421  ),
1422  map().pointMap()
1423  ).write();
1424  }
1425  }
1426 
1427  Info<< "End\n" << endl;
1428 
1429  return 0;
1430 }
1431 
1432 
1433 // ************************************************************************* //
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:51
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited (ie. from ordered back to original cell label)...
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:514
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:703
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:853
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1110
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: SortableList.H:56
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1072
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:41
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values (not the indices) of a list.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
const word dictName("faMeshDefinition")
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:402
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:888
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:639
static autoPtr< renumberMethod > New(const dictionary &renumberDict)
Return a reference to the selected renumbering method.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:365
void resetPrimitives(autoPtr< pointField > &&points, autoPtr< faceList > &&faces, autoPtr< labelList > &&owner, autoPtr< labelList > &&neighbour, const labelUList &patchSizes, const labelUList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:709
const cellList & cells() const
Field reading functions for post-processing utilities.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
label nFaces() const noexcept
Number of mesh faces.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
label nInternal1Edges() const
Internal edges using 0 or 1 boundary point.
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Abstract base class for renumbering.
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:847
void setSize(const label n)
Alias for resize()
Definition: List.H:289
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:964
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
virtual void resetAddressing(const labelUList &addr, const bool flipMapValue)
Reset addressing - use uniform flip map value.
Definition: faceZone.C:437
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:107
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:31
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:334
Cuthill-McKee renumbering.
label nInternalPoints() const noexcept
Points not on boundary.
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1104
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:376
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1293
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:977
label nInternalFaces() const noexcept
Number of internal faces.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1091
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:139
label nEdges() const
Number of mesh edges.
label nInternalEdges() const
Internal edges using 0,1 or 2 boundary points.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:646
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:75
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:29
static autoPtr< decompositionMethod > New(const dictionary &decompDict, const word &regionName="")
Return a reference to the selected decomposition method, optionally region-specific.
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
instantList times() const
Search the case for valid time directories.
Definition: TimePaths.C:142
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
const word & name() const
Return reference to name.
Definition: fvMesh.H:388
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:57
Direct mesh changes based on v1.3 polyTopoChange syntax.
const polyBoundaryMesh & patches
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:654
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:325
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
IOobject dictIO
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition: hexRef8Data.H:56
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:429
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
Foam::argList args(argc, argv)
ListType reorder(const labelUList &oldToNew, const ListType &input, const bool prune=false)
Reorder the elements of a list.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
virtual void updateMesh(const mapPolyMesh &morphMap)
Update any stored data for new labels. Not implemented.
Definition: topoSet.C:612
Foam::label startTime
label nInternal0Edges() const
Internal edges (i.e. not on boundary face) using no boundary point.
static void writeMaps(Ostream &os, const word &key, const labelListList &maps)
UIndirectList< bool > boolUIndList
UIndirectList of bools.
Definition: IndirectList.H:61
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:38
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157
Required Variables.