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 labelUList& elems
77 )
78 {
79  auto tfld = volScalarField::New
80  (
81  name,
82  mesh,
85  );
86  auto& fld = tfld.ref();
87 
88  forAll(fld, celli)
89  {
90  fld[celli] = elems[celli];
91  }
92 
93  return tfld;
94 }
95 
96 
97 // Calculate band of matrix
98 label getBand(const labelList& owner, const labelList& neighbour)
99 {
100  label band = 0;
101 
102  forAll(neighbour, facei)
103  {
104  label diff = neighbour[facei] - owner[facei];
105 
106  if (diff > band)
107  {
108  band = diff;
109  }
110  }
111  return band;
112 }
113 
114 
115 // Calculate band of matrix
116 void getBand
117 (
118  const bool calculateIntersect,
119  const label nCells,
120  const labelList& owner,
121  const labelList& neighbour,
122  label& bandwidth,
123  scalar& profile, // scalar to avoid overflow
124  scalar& sumSqrIntersect // scalar to avoid overflow
125 )
126 {
127  labelList cellBandwidth(nCells, Zero);
128  scalarField nIntersect(nCells, Zero);
129 
130  forAll(neighbour, facei)
131  {
132  label own = owner[facei];
133  label nei = neighbour[facei];
134 
135  // Note: mag not necessary for correct (upper-triangular) ordering.
136  label diff = nei-own;
137  cellBandwidth[nei] = max(cellBandwidth[nei], diff);
138  }
139 
140  bandwidth = max(cellBandwidth);
141 
142  // Do not use field algebra because of conversion label to scalar
143  profile = 0.0;
144  forAll(cellBandwidth, celli)
145  {
146  profile += 1.0*cellBandwidth[celli];
147  }
148 
149  sumSqrIntersect = 0.0;
150  if (calculateIntersect)
151  {
152  forAll(nIntersect, celli)
153  {
154  for (label colI = celli-cellBandwidth[celli]; colI <= celli; colI++)
155  {
156  nIntersect[colI] += 1.0;
157  }
158  }
159 
160  sumSqrIntersect = sum(Foam::sqr(nIntersect));
161  }
162 }
163 
164 
165 // Determine upper-triangular face order
166 labelList getFaceOrder
167 (
168  const primitiveMesh& mesh,
169  const labelList& cellOrder // New to old cell
170 )
171 {
172  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
173 
174  labelList oldToNewFace(mesh.nFaces(), -1);
175 
176  label newFacei = 0;
177 
178  labelList nbr;
179  labelList order;
180 
181  forAll(cellOrder, newCelli)
182  {
183  label oldCelli = cellOrder[newCelli];
184 
185  const cell& cFaces = mesh.cells()[oldCelli];
186 
187  // Neighbouring cells
188  nbr.setSize(cFaces.size());
189 
190  forAll(cFaces, i)
191  {
192  label facei = cFaces[i];
193 
194  if (mesh.isInternalFace(facei))
195  {
196  // Internal face. Get cell on other side.
197  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
198  if (nbrCelli == newCelli)
199  {
200  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
201  }
202 
203  if (newCelli < nbrCelli)
204  {
205  // Celli is master
206  nbr[i] = nbrCelli;
207  }
208  else
209  {
210  // nbrCell is master. Let it handle this face.
211  nbr[i] = -1;
212  }
213  }
214  else
215  {
216  // External face. Do later.
217  nbr[i] = -1;
218  }
219  }
220 
221  sortedOrder(nbr, order);
222 
223  for (const label index : order)
224  {
225  if (nbr[index] != -1)
226  {
227  oldToNewFace[cFaces[index]] = newFacei++;
228  }
229  }
230  }
231 
232  // Leave patch faces intact.
233  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
234  {
235  oldToNewFace[facei] = facei;
236  }
237 
238 
239  // Check done all faces.
240  forAll(oldToNewFace, facei)
241  {
242  if (oldToNewFace[facei] == -1)
243  {
245  << "Did not determine new position" << " for face " << facei
246  << abort(FatalError);
247  }
248  }
249 
250  return invert(mesh.nFaces(), oldToNewFace);
251 }
252 
253 
254 // Determine face order such that inside region faces are sorted
255 // upper-triangular but inbetween region faces are handled like boundary faces.
256 labelList getRegionFaceOrder
257 (
258  const primitiveMesh& mesh,
259  const labelList& cellOrder, // New to old cell
260  const labelList& cellToRegion // Old cell to region
261 )
262 {
263  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
264 
265  labelList oldToNewFace(mesh.nFaces(), -1);
266 
267  label newFacei = 0;
268 
269  label prevRegion = -1;
270 
271  forAll(cellOrder, newCelli)
272  {
273  label oldCelli = cellOrder[newCelli];
274 
275  if (cellToRegion[oldCelli] != prevRegion)
276  {
277  prevRegion = cellToRegion[oldCelli];
278  }
279 
280  const cell& cFaces = mesh.cells()[oldCelli];
281 
282  SortableList<label> nbr(cFaces.size());
283 
284  forAll(cFaces, i)
285  {
286  label facei = cFaces[i];
287 
288  if (mesh.isInternalFace(facei))
289  {
290  // Internal face. Get cell on other side.
291  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
292  if (nbrCelli == newCelli)
293  {
294  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
295  }
296 
297  if (cellToRegion[oldCelli] != cellToRegion[cellOrder[nbrCelli]])
298  {
299  // Treat like external face. Do later.
300  nbr[i] = -1;
301  }
302  else if (newCelli < nbrCelli)
303  {
304  // Celli is master
305  nbr[i] = nbrCelli;
306  }
307  else
308  {
309  // nbrCell is master. Let it handle this face.
310  nbr[i] = -1;
311  }
312  }
313  else
314  {
315  // External face. Do later.
316  nbr[i] = -1;
317  }
318  }
319 
320  nbr.sort();
321 
322  forAll(nbr, i)
323  {
324  if (nbr[i] != -1)
325  {
326  oldToNewFace[cFaces[nbr.indices()[i]]] = newFacei++;
327  }
328  }
329  }
330 
331  // Do region interfaces
332  label nRegions = max(cellToRegion)+1;
333  {
334  // Sort in increasing region
336 
337  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
338  {
339  label ownRegion = cellToRegion[mesh.faceOwner()[facei]];
340  label neiRegion = cellToRegion[mesh.faceNeighbour()[facei]];
341 
342  if (ownRegion != neiRegion)
343  {
344  sortKey[facei] =
345  min(ownRegion, neiRegion)*nRegions
346  +max(ownRegion, neiRegion);
347  }
348  }
349  sortKey.sort();
350 
351  // Extract.
352  label prevKey = -1;
353  forAll(sortKey, i)
354  {
355  label key = sortKey[i];
356 
357  if (key == labelMax)
358  {
359  break;
360  }
361 
362  if (prevKey != key)
363  {
364  prevKey = key;
365  }
366 
367  oldToNewFace[sortKey.indices()[i]] = newFacei++;
368  }
369  }
370 
371  // Leave patch faces intact.
372  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
373  {
374  oldToNewFace[facei] = facei;
375  }
376 
377 
378  // Check done all faces.
379  forAll(oldToNewFace, facei)
380  {
381  if (oldToNewFace[facei] == -1)
382  {
384  << "Did not determine new position"
385  << " for face " << facei
386  << abort(FatalError);
387  }
388  }
389 
390  return invert(mesh.nFaces(), oldToNewFace);
391 }
392 
393 
394 // cellOrder: old cell for every new cell
395 // faceOrder: old face for every new face. Ordering of boundary faces not
396 // changed.
397 autoPtr<mapPolyMesh> reorderMesh
398 (
399  polyMesh& mesh,
400  const labelList& cellOrder,
401  const labelList& faceOrder
402 )
403 {
404  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
405  labelList reverseFaceOrder(invert(faceOrder.size(), faceOrder));
406 
407  faceList newFaces(reorder(reverseFaceOrder, mesh.faces()));
408  labelList newOwner
409  (
410  renumber
411  (
412  reverseCellOrder,
413  reorder(reverseFaceOrder, mesh.faceOwner())
414  )
415  );
416  labelList newNeighbour
417  (
418  renumber
419  (
420  reverseCellOrder,
421  reorder(reverseFaceOrder, mesh.faceNeighbour())
422  )
423  );
424 
425  // Check if any faces need swapping.
426  labelHashSet flipFaceFlux(newOwner.size());
427  forAll(newNeighbour, facei)
428  {
429  label own = newOwner[facei];
430  label nei = newNeighbour[facei];
431 
432  if (nei < own)
433  {
434  newFaces[facei].flip();
435  std::swap(newOwner[facei], newNeighbour[facei]);
436  flipFaceFlux.insert(facei);
437  }
438  }
439 
441  labelList patchSizes(patches.size());
442  labelList patchStarts(patches.size());
443  labelList oldPatchNMeshPoints(patches.size());
444  labelListList patchPointMap(patches.size());
445 
446  forAll(patches, patchi)
447  {
448  patchSizes[patchi] = patches[patchi].size();
449  patchStarts[patchi] = patches[patchi].start();
450  oldPatchNMeshPoints[patchi] = patches[patchi].nPoints();
451  patchPointMap[patchi] = identity(patches[patchi].nPoints());
452  }
453 
455  (
456  autoPtr<pointField>(), // <- null: leaves points untouched
457  autoPtr<faceList>::New(std::move(newFaces)),
458  autoPtr<labelList>::New(std::move(newOwner)),
459  autoPtr<labelList>::New(std::move(newNeighbour)),
460  patchSizes,
461  patchStarts,
462  true
463  );
464 
465 
466  // Re-do the faceZones
467  {
468  faceZoneMesh& faceZones = mesh.faceZones();
469  faceZones.clearAddressing();
470  forAll(faceZones, zoneI)
471  {
472  faceZone& fZone = faceZones[zoneI];
473  labelList newAddressing(fZone.size());
474  boolList newFlipMap(fZone.size());
475  forAll(fZone, i)
476  {
477  label oldFacei = fZone[i];
478  newAddressing[i] = reverseFaceOrder[oldFacei];
479  if (flipFaceFlux.found(newAddressing[i]))
480  {
481  newFlipMap[i] = !fZone.flipMap()[i];
482  }
483  else
484  {
485  newFlipMap[i] = fZone.flipMap()[i];
486  }
487  }
488  labelList newToOld(sortedOrder(newAddressing));
489  fZone.resetAddressing
490  (
491  labelUIndList(newAddressing, newToOld)(),
492  boolUIndList(newFlipMap, newToOld)()
493  );
494  }
495  }
496  // Re-do the cellZones
497  {
498  cellZoneMesh& cellZones = mesh.cellZones();
499  cellZones.clearAddressing();
500  forAll(cellZones, zoneI)
501  {
502  cellZones[zoneI] = labelUIndList
503  (
504  reverseCellOrder,
505  cellZones[zoneI]
506  )();
507  Foam::sort(cellZones[zoneI]);
508  }
509  }
510 
511 
513  (
514  mesh, // const polyMesh& mesh,
515  mesh.nPoints(), // nOldPoints,
516  mesh.nFaces(), // nOldFaces,
517  mesh.nCells(), // nOldCells,
518  identity(mesh.nPoints()), // pointMap,
519  List<objectMap>(), // pointsFromPoints,
520  faceOrder, // faceMap,
521  List<objectMap>(), // facesFromPoints,
522  List<objectMap>(), // facesFromEdges,
523  List<objectMap>(), // facesFromFaces,
524  cellOrder, // cellMap,
525  List<objectMap>(), // cellsFromPoints,
526  List<objectMap>(), // cellsFromEdges,
527  List<objectMap>(), // cellsFromFaces,
528  List<objectMap>(), // cellsFromCells,
529  identity(mesh.nPoints()), // reversePointMap,
530  reverseFaceOrder, // reverseFaceMap,
531  reverseCellOrder, // reverseCellMap,
532  flipFaceFlux, // flipFaceFlux,
533  patchPointMap, // patchPointMap,
534  labelListList(), // pointZoneMap,
535  labelListList(), // faceZonePointMap,
536  labelListList(), // faceZoneFaceMap,
537  labelListList(), // cellZoneMap,
538  pointField(), // preMotionPoints,
539  patchStarts, // oldPatchStarts,
540  oldPatchNMeshPoints, // oldPatchNMeshPoints
541  autoPtr<scalarField>() // oldCellVolumes
542  );
543 }
544 
545 
546 // Return new to old cell numbering
547 labelList regionRenumber
548 (
549  const renumberMethod& method,
550  const fvMesh& mesh,
551  const labelList& cellToRegion
552 )
553 {
554  Info<< "Determining cell order:" << endl;
555 
556  labelList cellOrder(cellToRegion.size());
557 
558  label nRegions = max(cellToRegion)+1;
559 
560  labelListList regionToCells(invertOneToMany(nRegions, cellToRegion));
561 
562  label celli = 0;
563 
564  forAll(regionToCells, regioni)
565  {
566  Info<< " region " << regioni << " starts at " << celli << endl;
567 
568  // Make sure no parallel comms
569  const bool oldParRun = UPstream::parRun(false);
570 
571  // Per region do a reordering.
572  fvMeshSubset subsetter(mesh, regioni, cellToRegion);
573 
574  const fvMesh& subMesh = subsetter.subMesh();
575 
576  labelList subCellOrder = method.renumber
577  (
578  subMesh,
579  subMesh.cellCentres()
580  );
581 
582  UPstream::parRun(oldParRun); // Restore parallel state
583 
584 
585  const labelList& cellMap = subsetter.cellMap();
586 
587  forAll(subCellOrder, i)
588  {
589  cellOrder[celli++] = cellMap[subCellOrder[i]];
590  }
591  }
592  Info<< endl;
593 
594  return cellOrder;
595 }
596 
597 
598 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
599 
600 int main(int argc, char *argv[])
601 {
603  (
604  "Renumber mesh cells to reduce the bandwidth"
605  );
606 
607  #include "addAllRegionOptions.H"
608  #include "addOverwriteOption.H"
609  #include "addTimeOptions.H"
610 
611  argList::addOption("dict", "file", "Alternative renumberMeshDict");
612 
614  (
615  "frontWidth",
616  "Calculate the rms of the front-width"
617  );
618 
619  argList::noFunctionObjects(); // Never use function objects
620 
621  #include "setRootCase.H"
622  #include "createTime.H"
623  #include "getAllRegionOptions.H"
624 
625  // Force linker to include zoltan symbols. This section is only needed since
626  // Zoltan is a static library
627  #ifdef HAVE_ZOLTAN
628  Info<< "renumberMesh built with zoltan support." << nl << endl;
629  (void)zoltanRenumber::typeName;
630  #endif
631 
632  // Get times list
633  instantList Times = runTime.times();
634 
635  // Set startTime and endTime depending on -time and -latestTime options
636  #include "checkTimeOptions.H"
637 
639 
640  #include "createNamedMeshes.H"
641 
642  const bool readDict = args.found("dict");
643  const bool doFrontWidth = args.found("frontWidth");
644  const bool overwrite = args.found("overwrite");
645 
646 
647  for (fvMesh& mesh : meshes)
648  {
649  // Reset time in case of multiple meshes and not overwrite
651 
652  const word oldInstance = mesh.pointsInstance();
653 
654  label band;
655  scalar profile;
656  scalar sumSqrIntersect;
657  getBand
658  (
659  doFrontWidth,
660  mesh.nCells(),
661  mesh.faceOwner(),
663  band,
664  profile,
665  sumSqrIntersect
666  );
667 
668  reduce(band, maxOp<label>());
669  reduce(profile, sumOp<scalar>());
670  scalar rmsFrontwidth = Foam::sqrt
671  (
673  (
674  sumSqrIntersect,
675  sumOp<scalar>()
677  );
678 
679  Info<< "Mesh " << mesh.name()
680  << " size: " << mesh.globalData().nTotalCells() << nl
681  << "Before renumbering :" << nl
682  << " band : " << band << nl
683  << " profile : " << profile << nl;
684 
685  if (doFrontWidth)
686  {
687  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
688  }
689 
690  Info<< endl;
691 
692  bool sortCoupledFaceCells = false;
693  bool writeMaps = false;
694  bool orderPoints = false;
695  label blockSize = 0;
696 
697  // Construct renumberMethod
698  autoPtr<IOdictionary> renumberDictPtr;
699  autoPtr<renumberMethod> renumberPtr;
700 
701  if (readDict)
702  {
703  const word dictName("renumberMeshDict");
704  #include "setSystemMeshDictionaryIO.H"
705 
706  Info<< "Renumber according to " << dictIO.name() << nl << endl;
707 
708  renumberDictPtr.reset(new IOdictionary(dictIO));
709  const IOdictionary& renumberDict = renumberDictPtr();
710 
711  renumberPtr = renumberMethod::New(renumberDict);
712 
713  sortCoupledFaceCells = renumberDict.getOrDefault
714  (
715  "sortCoupledFaceCells",
716  false
717  );
718  if (sortCoupledFaceCells)
719  {
720  Info<< "Sorting cells on coupled boundaries to be last." << nl
721  << endl;
722  }
723 
724  blockSize = renumberDict.getOrDefault("blockSize", 0);
725  if (blockSize > 0)
726  {
727  Info<< "Ordering cells into regions of size " << blockSize
728  << " (using decomposition);"
729  << " ordering faces into region-internal"
730  << " and region-external."
731  << nl << endl;
732 
733  if (blockSize < 0 || blockSize >= mesh.nCells())
734  {
736  << "Block size " << blockSize
737  << " should be positive integer"
738  << " and less than the number of cells in the mesh."
739  << exit(FatalError);
740  }
741  }
742 
743  orderPoints = renumberDict.getOrDefault("orderPoints", false);
744  if (orderPoints)
745  {
746  Info<< "Ordering points into internal and boundary points."
747  << nl
748  << endl;
749  }
750 
751  renumberDict.readEntry("writeMaps", writeMaps);
752  if (writeMaps)
753  {
754  Info<< "Writing renumber maps (new to old) to polyMesh." << nl
755  << endl;
756  }
757  }
758  else
759  {
760  Info<< "Using default renumberMethod." << nl << endl;
761  dictionary renumberDict;
762  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
763  }
764 
765  Info<< "Selecting renumberMethod " << renumberPtr().type() << nl
766  << endl;
767 
768 
769 
770  // Read parallel reconstruct maps
771  labelIOList cellProcAddressing
772  (
773  IOobject
774  (
775  "cellProcAddressing",
778  mesh,
781  ),
782  labelList(0)
783  );
784 
785  labelIOList faceProcAddressing
786  (
787  IOobject
788  (
789  "faceProcAddressing",
792  mesh,
795  ),
796  labelList(0)
797  );
798  labelIOList pointProcAddressing
799  (
800  IOobject
801  (
802  "pointProcAddressing",
805  mesh,
808  ),
809  labelList(0)
810  );
811  labelIOList boundaryProcAddressing
812  (
813  IOobject
814  (
815  "boundaryProcAddressing",
818  mesh,
821  ),
822  labelList(0)
823  );
824 
825 
826  // Read objects in time directory
827  IOobjectList objects(mesh, runTime.timeName());
828 
829 
830  // Read vol fields.
831 
833  ReadFields(mesh, objects, vsFlds);
834 
836  ReadFields(mesh, objects, vvFlds);
837 
839  ReadFields(mesh, objects, vstFlds);
840 
841  PtrList<volSymmTensorField> vsymtFlds;
842  ReadFields(mesh, objects, vsymtFlds);
843 
845  ReadFields(mesh, objects, vtFlds);
846 
847 
848  // Read surface fields.
849 
851  ReadFields(mesh, objects, ssFlds);
852 
854  ReadFields(mesh, objects, svFlds);
855 
857  ReadFields(mesh, objects, sstFlds);
858 
860  ReadFields(mesh, objects, ssymtFlds);
861 
863  ReadFields(mesh, objects, stFlds);
864 
865 
866  // Read point fields.
867 
869  ReadFields(pointMesh::New(mesh), objects, psFlds);
870 
872  ReadFields(pointMesh::New(mesh), objects, pvFlds);
873 
875  ReadFields(pointMesh::New(mesh), objects, pstFlds);
876 
878  ReadFields(pointMesh::New(mesh), objects, psymtFlds);
879 
881  ReadFields(pointMesh::New(mesh), objects, ptFlds);
882 
883 
884  // Read sets
885  PtrList<cellSet> cellSets;
886  PtrList<faceSet> faceSets;
887  PtrList<pointSet> pointSets;
888  {
889  // Read sets
890  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
891  ReadFields(objects, cellSets);
892  ReadFields(objects, faceSets);
893  ReadFields(objects, pointSets);
894  }
895 
896 
897  Info<< endl;
898 
899  // From renumbering:
900  // - from new cell/face back to original cell/face
901  labelList cellOrder;
902  labelList faceOrder;
903 
904  if (blockSize > 0)
905  {
906  // Renumbering in two phases. Should be done in one so mapping of
907  // fields is done correctly!
908 
909  label nBlocks = mesh.nCells()/blockSize;
910  Info<< "nBlocks = " << nBlocks << endl;
911 
912  // Read decompositionMethod dictionary
913  dictionary decomposeDict(renumberDictPtr().subDict("blockCoeffs"));
914  decomposeDict.set("numberOfSubdomains", nBlocks);
915 
916  const bool oldParRun = UPstream::parRun(false);
917 
919  (
920  decomposeDict
921  );
922 
923  labelList cellToRegion
924  (
925  decomposePtr().decompose
926  (
927  mesh,
928  mesh.cellCentres()
929  )
930  );
931 
932  UPstream::parRun(oldParRun); // Restore parallel state
933 
934 
935  // For debugging: write out region
936  createScalarField
937  (
938  mesh,
939  "cellDist",
940  cellToRegion
941  )().write();
942 
943  Info<< nl << "Written decomposition as volScalarField to "
944  << "cellDist for use in postprocessing."
945  << nl << endl;
946 
947 
948  cellOrder = regionRenumber(renumberPtr(), mesh, cellToRegion);
949 
950  // Determine new to old face order with new cell numbering
951  faceOrder = getRegionFaceOrder
952  (
953  mesh,
954  cellOrder,
955  cellToRegion
956  );
957  }
958  else
959  {
960  // Determines sorted back to original cell ordering
961  cellOrder = renumberPtr().renumber
962  (
963  mesh,
964  mesh.cellCentres()
965  );
966 
967  if (sortCoupledFaceCells)
968  {
969  // Change order so all coupled patch faceCells are at the end.
971 
972  // Collect all boundary cells on coupled patches
973  label nBndCells = 0;
974  forAll(pbm, patchi)
975  {
976  if (pbm[patchi].coupled())
977  {
978  nBndCells += pbm[patchi].size();
979  }
980  }
981 
982  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
983 
984  labelList bndCells(nBndCells);
985  labelList bndCellMap(nBndCells);
986  nBndCells = 0;
987  forAll(pbm, patchi)
988  {
989  if (pbm[patchi].coupled())
990  {
991  const labelUList& faceCells = pbm[patchi].faceCells();
992  forAll(faceCells, i)
993  {
994  label celli = faceCells[i];
995 
996  if (reverseCellOrder[celli] != -1)
997  {
998  bndCells[nBndCells] = celli;
999  bndCellMap[nBndCells++] =
1000  reverseCellOrder[celli];
1001  reverseCellOrder[celli] = -1;
1002  }
1003  }
1004  }
1005  }
1006  bndCells.setSize(nBndCells);
1007  bndCellMap.setSize(nBndCells);
1008 
1009  // Sort
1010  labelList order(sortedOrder(bndCellMap));
1011 
1012  // Redo newReverseCellOrder
1013  labelList newReverseCellOrder(mesh.nCells(), -1);
1014 
1015  label sortedI = mesh.nCells();
1016  forAllReverse(order, i)
1017  {
1018  label origCelli = bndCells[order[i]];
1019  newReverseCellOrder[origCelli] = --sortedI;
1020  }
1021 
1022  Info<< "Ordered all " << nBndCells
1023  << " cells with a coupled face"
1024  << " to the end of the cell list, starting at " << sortedI
1025  << endl;
1026 
1027  // Compact
1028  sortedI = 0;
1029  forAll(cellOrder, newCelli)
1030  {
1031  label origCelli = cellOrder[newCelli];
1032  if (newReverseCellOrder[origCelli] == -1)
1033  {
1034  newReverseCellOrder[origCelli] = sortedI++;
1035  }
1036  }
1037 
1038  // Update sorted back to original (unsorted) map
1039  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
1040  }
1041 
1042 
1043  // Determine new to old face order with new cell numbering
1044  faceOrder = getFaceOrder
1045  (
1046  mesh,
1047  cellOrder // New to old cell
1048  );
1049  }
1050 
1051 
1052  if (!overwrite)
1053  {
1054  ++runTime;
1055  }
1056 
1057 
1058  // Change the mesh.
1059  autoPtr<mapPolyMesh> map = reorderMesh(mesh, cellOrder, faceOrder);
1060 
1061 
1062  if (orderPoints)
1063  {
1064  polyTopoChange meshMod(mesh);
1065  autoPtr<mapPolyMesh> pointOrderMap = meshMod.changeMesh
1066  (
1067  mesh,
1068  false, // inflate
1069  true, // syncParallel
1070  false, // orderCells
1071  orderPoints // orderPoints
1072  );
1073 
1074  // Combine point reordering into map.
1075  const_cast<labelList&>(map().pointMap()) = labelUIndList
1076  (
1077  map().pointMap(),
1078  pointOrderMap().pointMap()
1079  )();
1080 
1082  (
1083  pointOrderMap().reversePointMap(),
1084  const_cast<labelList&>(map().reversePointMap())
1085  );
1086  }
1087 
1088 
1089  // Update fields
1090  mesh.updateMesh(map());
1091 
1092  // Update proc maps
1093  if (cellProcAddressing.headerOk())
1094  {
1095  if (returnReduceAnd(cellProcAddressing.size() == mesh.nCells()))
1096  {
1097  Info<< "Renumbering processor cell decomposition map "
1098  << cellProcAddressing.name() << endl;
1099 
1100  cellProcAddressing = labelList
1101  (
1102  labelUIndList(cellProcAddressing, map().cellMap())
1103  );
1104  }
1105  else
1106  {
1107  Info<< "Not writing inconsistent processor cell decomposition"
1108  << " map " << cellProcAddressing.filePath() << endl;
1109  cellProcAddressing.writeOpt(IOobject::NO_WRITE);
1110  }
1111  }
1112  else
1113  {
1114  cellProcAddressing.writeOpt(IOobject::NO_WRITE);
1115  }
1116 
1117  if (faceProcAddressing.headerOk())
1118  {
1119  if (returnReduceAnd(faceProcAddressing.size() == mesh.nFaces()))
1120  {
1121  Info<< "Renumbering processor face decomposition map "
1122  << faceProcAddressing.name() << endl;
1123 
1124  faceProcAddressing = labelList
1125  (
1126  labelUIndList(faceProcAddressing, map().faceMap())
1127  );
1128 
1129  // Detect any flips.
1130  const labelHashSet& fff = map().flipFaceFlux();
1131  for (const label facei : fff)
1132  {
1133  label masterFacei = faceProcAddressing[facei];
1134 
1135  faceProcAddressing[facei] = -masterFacei;
1136 
1137  if (masterFacei == 0)
1138  {
1140  << " masterFacei:" << masterFacei
1141  << exit(FatalError);
1142  }
1143  }
1144  }
1145  else
1146  {
1147  Info<< "Not writing inconsistent processor face decomposition"
1148  << " map " << faceProcAddressing.filePath() << endl;
1149  faceProcAddressing.writeOpt(IOobject::NO_WRITE);
1150  }
1151  }
1152  else
1153  {
1154  faceProcAddressing.writeOpt(IOobject::NO_WRITE);
1155  }
1156 
1157  if (pointProcAddressing.headerOk())
1158  {
1159  if (returnReduceAnd(pointProcAddressing.size() == mesh.nPoints()))
1160  {
1161  Info<< "Renumbering processor point decomposition map "
1162  << pointProcAddressing.name() << endl;
1163 
1164  pointProcAddressing = labelList
1165  (
1166  labelUIndList(pointProcAddressing, map().pointMap())
1167  );
1168  }
1169  else
1170  {
1171  Info<< "Not writing inconsistent processor point decomposition"
1172  << " map " << pointProcAddressing.filePath() << endl;
1173  pointProcAddressing.writeOpt(IOobject::NO_WRITE);
1174  }
1175  }
1176  else
1177  {
1178  pointProcAddressing.writeOpt(IOobject::NO_WRITE);
1179  }
1180 
1181  if (boundaryProcAddressing.headerOk())
1182  {
1183  if
1184  (
1186  (
1187  boundaryProcAddressing.size() == mesh.boundaryMesh().size()
1188  )
1189  )
1190  {
1191  // No renumbering needed
1192  }
1193  else
1194  {
1195  Info<< "Not writing inconsistent processor patch decomposition"
1196  << " map " << boundaryProcAddressing.filePath() << endl;
1197  boundaryProcAddressing.writeOpt(IOobject::NO_WRITE);
1198  }
1199  }
1200  else
1201  {
1202  boundaryProcAddressing.writeOpt(IOobject::NO_WRITE);
1203  }
1204 
1205 
1206 
1207 
1208  // Move mesh (since morphing might not do this)
1209  if (map().hasMotionPoints())
1210  {
1211  mesh.movePoints(map().preMotionPoints());
1212  }
1213 
1214 
1215  {
1216  label band;
1217  scalar profile;
1218  scalar sumSqrIntersect;
1219  getBand
1220  (
1221  doFrontWidth,
1222  mesh.nCells(),
1223  mesh.faceOwner(),
1224  mesh.faceNeighbour(),
1225  band,
1226  profile,
1227  sumSqrIntersect
1228  );
1229  reduce(band, maxOp<label>());
1230  reduce(profile, sumOp<scalar>());
1231  scalar rmsFrontwidth = Foam::sqrt
1232  (
1233  returnReduce
1234  (
1235  sumSqrIntersect,
1236  sumOp<scalar>()
1238  );
1239 
1240  Info<< "After renumbering :" << nl
1241  << " band : " << band << nl
1242  << " profile : " << profile << nl;
1243 
1244  if (doFrontWidth)
1245  {
1246 
1247  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1248  }
1249 
1250  Info<< endl;
1251  }
1252 
1253  if (orderPoints)
1254  {
1255  // Force edge calculation (since only reason that points would
1256  // need to be sorted)
1257  (void)mesh.edges();
1258 
1259  label nTotPoints = returnReduce
1260  (
1261  mesh.nPoints(),
1262  sumOp<label>()
1263  );
1264  label nTotIntPoints = returnReduce
1265  (
1267  sumOp<label>()
1268  );
1269 
1270  label nTotEdges = returnReduce
1271  (
1272  mesh.nEdges(),
1273  sumOp<label>()
1274  );
1275  label nTotIntEdges = returnReduce
1276  (
1277  mesh.nInternalEdges(),
1278  sumOp<label>()
1279  );
1280  label nTotInt0Edges = returnReduce
1281  (
1283  sumOp<label>()
1284  );
1285  label nTotInt1Edges = returnReduce
1286  (
1288  sumOp<label>()
1289  );
1290 
1291  Info<< "Points:" << nl
1292  << " total : " << nTotPoints << nl
1293  << " internal: " << nTotIntPoints << nl
1294  << " boundary: " << nTotPoints-nTotIntPoints << nl
1295  << "Edges:" << nl
1296  << " total : " << nTotEdges << nl
1297  << " internal: " << nTotIntEdges << nl
1298  << " internal using 0 boundary points: "
1299  << nTotInt0Edges << nl
1300  << " internal using 1 boundary points: "
1301  << nTotInt1Edges-nTotInt0Edges << nl
1302  << " internal using 2 boundary points: "
1303  << nTotIntEdges-nTotInt1Edges << nl
1304  << " boundary: " << nTotEdges-nTotIntEdges << nl
1305  << endl;
1306  }
1307 
1308 
1309  if (overwrite)
1310  {
1311  mesh.setInstance(oldInstance);
1312  }
1313  else
1314  {
1316  }
1317 
1318 
1319  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1320 
1321  // Remove old procAddressing files
1323 
1324  // Update refinement data
1325  hexRef8Data refData
1326  (
1327  IOobject
1328  (
1329  "dummy",
1330  mesh.facesInstance(),
1332  mesh,
1336  )
1337  );
1338  refData.updateMesh(map());
1339  refData.write();
1340 
1341  // Update sets
1342  topoSet::updateMesh(mesh.facesInstance(), map(), cellSets);
1343  topoSet::updateMesh(mesh.facesInstance(), map(), faceSets);
1344  topoSet::updateMesh(mesh.facesInstance(), map(), pointSets);
1345 
1346  mesh.write();
1347 
1348  if (writeMaps)
1349  {
1350  // For debugging: write out region
1351  createScalarField
1352  (
1353  mesh,
1354  "origCellID",
1355  map().cellMap()
1356  )().write();
1357 
1358  createScalarField
1359  (
1360  mesh,
1361  "cellID",
1362  identity(mesh.nCells())
1363  )().write();
1364 
1365  Info<< nl
1366  << "Written current cellID and origCellID as volScalarField"
1367  << " for use in postprocessing." << nl << endl;
1368 
1369  labelIOList
1370  (
1371  IOobject
1372  (
1373  "cellMap",
1374  mesh.facesInstance(),
1376  mesh,
1380  ),
1381  map().cellMap()
1382  ).write();
1383 
1384  labelIOList
1385  (
1386  IOobject
1387  (
1388  "faceMap",
1389  mesh.facesInstance(),
1391  mesh,
1395  ),
1396  map().faceMap()
1397  ).write();
1398 
1399  labelIOList
1400  (
1401  IOobject
1402  (
1403  "pointMap",
1404  mesh.facesInstance(),
1406  mesh,
1410  ),
1411  map().pointMap()
1412  ).write();
1413  }
1414  }
1415 
1416  Info<< "End\n" << endl;
1417 
1418  return 0;
1419 }
1420 
1421 
1422 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
Foam::surfaceFields.
const polyBoundaryMesh & pbm
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:547
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:815
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:859
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:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
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:195
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: SortableList.H:56
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. Registered with typeName.
Definition: MeshObject.C:53
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:50
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
const word dictName("faMeshDefinition")
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:918
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
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:374
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:704
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
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.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Required Classes.
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)
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:421
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:853
IOList< label > labelIOList
IO for a List of label.
Definition: labelIOList.H:32
void setSize(const label n)
Alias for resize()
Definition: List.H:316
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:994
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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), works like std::iota() but returning a...
Definition: labelLists.C:44
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
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:296
Cuthill-McKee renumbering.
label nInternalPoints() const noexcept
Points not on boundary.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
Reading is optional [identical to LAZY_READ].
virtual void resetAddressing(faceZone &&zn)
Move reset addressing and flip map from another zone.
Definition: faceZone.C:477
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:385
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:935
label nInternalFaces() const noexcept
Number of internal faces.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
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.
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1102
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 tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
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.
instantList times() const
Search the case for valid time directories.
Definition: TimePaths.C:118
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:387
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:78
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
Direct mesh changes based on v1.3 polyTopoChange syntax.
const polyBoundaryMesh & patches
Nothing to be read.
Required Classes.
Automatically write from objectRegistry::writeObject()
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
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:367
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:74
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:437
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
bool coupled
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:172
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)
Do not request registration (bool: false)
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
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
Required Variables.