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-2025 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 (Cuthill-McKee) or the method
38  specified by the -renumber-method option, but will read
39  system/renumberMeshDict if -dict option is present
40 
41 Usage
42  \b renumberMesh [OPTIONS]
43 
44  Options:
45  - \par -allRegions
46  Use all regions in regionProperties
47 
48  - \par -case <dir>
49  Specify case directory to use (instead of the cwd).
50 
51  - \par -constant
52  Include the 'constant/' dir in the times list.
53 
54  - \par -decompose
55  Aggregate initially with a decomposition method (serial only)
56 
57  - \par -decomposeParDict <file>
58  Use specified file for decomposePar dictionary.
59 
60  - \par -dict <file>
61  Use specified file for renumberMeshDict dictionary.
62 
63  - \par -dry-run
64  Test only
65 
66  - \par -frontWidth
67  Calculate the rms of the front-width
68 
69  - \par -latestTime
70  Select the latest time.
71 
72  - \par -lib <name>
73  Additional library or library list to load (can be used multiple times).
74 
75  - \par -no-fields
76  Suppress renumber of fields
77 
78  - \par -noZero
79  Exclude the \a 0 dir from the times list.
80 
81  - \par -overwrite
82  Overwrite existing mesh/results files
83 
84  - \par -parallel
85  Run in parallel
86 
87  - \par -region <regionName>
88  Renumber named region.
89 
90  - \par -regions <wordRes>
91  Renumber named regions.
92 
93  - \par -renumber-coeffs <string-content>
94  String to create renumber dictionary contents.
95 
96  - \par -renumber-method <name>
97  Specify renumber method (default: CuthillMcKee) without dictionary
98 
99  - \par -time <value>
100  Specify time to select
101 
102  - \par -verbose
103  Additional verbosity.
104 
105  - \par -doc
106  Display documentation in browser.
107 
108  - \par -doc-source
109  Display source code in browser.
110 
111  - \par -help
112  Display short help and exit.
113 
114  - \par -help-man
115  Display full help (manpage format) and exit.
116 
117  - \par -help-notes
118  Display help notes (description) and exit.
119 
120  - \par -help-full
121  Display full help and exit.
122 
123 \*---------------------------------------------------------------------------*/
124 
125 #include "argList.H"
126 #include "clockTime.H"
127 #include "timeSelector.H"
128 #include "IOobjectList.H"
129 #include "fvMesh.H"
130 #include "polyTopoChange.H"
131 #include "ReadFields.H"
132 #include "volFields.H"
133 #include "surfaceFields.H"
134 #include "SortableList.H"
135 #include "decompositionMethod.H"
136 #include "decompositionModel.H"
137 #include "renumberMethod.H"
139 #include "CuthillMcKeeRenumber.H"
140 #include "fvMeshSubset.H"
141 #include "cellSet.H"
142 #include "faceSet.H"
143 #include "pointSet.H"
144 #include "processorMeshes.H"
145 #include "hexRef8Data.H"
146 #include "regionProperties.H"
147 #include "polyMeshTools.H"
148 #include "subsetAdjacency.H"
149 
150 using namespace Foam;
151 
152 // Slightly messy way of handling timing, but since the timing points
153 // are scattered between 'main()' and other local functions...
154 
156 
157 // Timing categories
158 enum TimingType
159 {
160  READ_MESH, // Reading mesh
161  READ_FIELDS, // Reading fields
162  DECOMPOSE, // Domain decomposition (if any)
163  CELL_CELLS, // globalMeshData::calcCellCells
164  RENUMBER, // The renumberMethod
165  REORDER, // Mesh reordering (topoChange)
166  WRITING, // Writing mesh/fields
167 };
168 FixedList<double, 8> timings;
169 
170 
171 // Create named field from labelList for post-processing
172 tmp<volScalarField> createScalarField
173 (
174  const fvMesh& mesh,
175  const word& name,
176  const labelUList& elems
177 )
178 {
179  auto tfld = volScalarField::New
180  (
181  name,
182  mesh,
185  );
186  auto& fld = tfld.ref();
187 
188  forAll(elems, celli)
189  {
190  fld[celli] = elems[celli];
191  }
192 
193  fld.correctBoundaryConditions();
194  return tfld;
195 }
196 
197 
198 // Calculate band of mesh
199 // label getBand(const labelUList& owner, const labelUList& neighbour)
200 // {
201 // label bandwidth = 0;
202 //
203 // forAll(neighbour, facei)
204 // {
205 // const label width = neighbour[facei] - owner[facei];
206 //
207 // if (bandwidth < width)
208 // {
209 // bandwidth = width;
210 // }
211 // }
212 // return bandwidth;
213 // }
214 
215 
216 // Calculate band and profile of matrix. Profile is scalar to avoid overflow
217 Tuple2<label, scalar> getBand
218 (
219  const CompactListList<label>& mat
220 )
221 {
222  Tuple2<label, scalar> metrics(0, 0);
223 
224  auto& bandwidth = metrics.first();
225  auto& profile = metrics.second();
226 
227  forAll(mat, celli)
228  {
229  const auto& neighbours = mat[celli];
230 
231  const label nNbr = neighbours.size();
232 
233  if (nNbr)
234  {
235  // Max distance
236  const label width = (neighbours[nNbr-1] - celli);
237 
238  if (bandwidth < width)
239  {
240  bandwidth = width;
241  }
242 
243  profile += scalar(width);
244  }
245  }
246 
247  return metrics;
248 }
249 
250 
251 // Calculate band of matrix
252 void getBand
253 (
254  const bool calculateIntersect,
255  const label nCells,
256  const labelUList& owner,
257  const labelUList& neighbour,
258  label& bandwidth,
259  scalar& profile, // scalar to avoid overflow
260  scalar& sumSqrIntersect // scalar to avoid overflow
261 )
262 {
263  labelList cellBandwidth(nCells, Foam::zero{});
264 
265  bandwidth = 0;
266 
267  forAll(neighbour, facei)
268  {
269  const label own = owner[facei];
270  const label nei = neighbour[facei];
271 
272  // Note: mag not necessary for correct (upper-triangular) ordering.
273  const label width = nei - own;
274 
275  if (cellBandwidth[nei] < width)
276  {
277  cellBandwidth[nei] = width;
278 
279  if (bandwidth < width)
280  {
281  bandwidth = width;
282  }
283  }
284  }
285 
286  // Do not use field algebra because of conversion label to scalar
287  profile = 0;
288  for (const label width : cellBandwidth)
289  {
290  profile += scalar(width);
291  }
292 
293  sumSqrIntersect = 0;
294  if (calculateIntersect)
295  {
296  scalarField nIntersect(nCells, Foam::zero{});
297 
298  forAll(nIntersect, celli)
299  {
300  for (label colI = celli-cellBandwidth[celli]; colI <= celli; colI++)
301  {
302  nIntersect[colI] += scalar(1);
303  }
304  }
305 
306  sumSqrIntersect = sum(Foam::sqr(nIntersect));
307  }
308 }
309 
310 
311 // Determine upper-triangular face order
312 labelList getFaceOrder
313 (
314  const primitiveMesh& mesh,
315  const labelList& cellOrder // New to old cell
316 )
317 {
318  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
319 
320  labelList oldToNewFace(mesh.nFaces(), -1);
321 
322  label newFacei = 0;
323 
324  DynamicList<label> nbr(64);
325  DynamicList<label> order(64);
326 
327  forAll(cellOrder, newCelli)
328  {
329  const label oldCelli = cellOrder[newCelli];
330 
331  const cell& cFaces = mesh.cells()[oldCelli];
332 
333  // Neighbouring cells
334  nbr.clear();
335 
336  for (const label facei : cFaces)
337  {
338  label nbrCelli = -1;
339 
340  if (mesh.isInternalFace(facei))
341  {
342  // Internal face. Get cell on other side.
343  nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
344  if (nbrCelli == newCelli)
345  {
346  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
347  }
348 
349  // The nbrCell is actually the master (let it handle the face)
350  if (nbrCelli <= newCelli)
351  {
352  nbrCelli = -1;
353  }
354  }
355 
356  nbr.push_back(nbrCelli);
357  }
358 
359  Foam::sortedOrder(nbr, order);
360 
361  for (const label index : order)
362  {
363  if (nbr[index] >= 0)
364  {
365  oldToNewFace[cFaces[index]] = newFacei++;
366  }
367  }
368  }
369 
370  // Leave patch faces intact.
371  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
372  {
373  oldToNewFace[facei] = facei;
374  }
375 
376 
377  // Check done all faces.
378  forAll(oldToNewFace, facei)
379  {
380  if (oldToNewFace[facei] == -1)
381  {
383  << "Did not determine new position" << " for face " << facei
384  << abort(FatalError);
385  }
386  }
387 
388  return invert(mesh.nFaces(), oldToNewFace);
389 }
390 
391 
392 // Determine face order such that inside region faces are sorted
393 // upper-triangular but inbetween region faces are handled like boundary faces.
394 labelList getRegionFaceOrder
395 (
396  const primitiveMesh& mesh,
397  const labelList& cellOrder, // New to old cell
398  const labelList& cellToRegion // Old cell to region
399 )
400 {
401  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
402 
403  labelList oldToNewFace(mesh.nFaces(), -1);
404 
405  label newFacei = 0;
406 
407  DynamicList<label> nbr(64);
408  DynamicList<label> order(64);
409 
410  forAll(cellOrder, newCelli)
411  {
412  const label oldCelli = cellOrder[newCelli];
413 
414  const cell& cFaces = mesh.cells()[oldCelli];
415 
416  // Neighbouring cells
417  nbr.clear();
418 
419  for (const label facei : cFaces)
420  {
421  label nbrCelli = -1;
422 
423  if (mesh.isInternalFace(facei))
424  {
425  // Internal face. Get cell on other side.
426  nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
427  if (nbrCelli == newCelli)
428  {
429  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
430  }
431 
432  // The nbrCell is actually the master (let it handle the face)
433  if (nbrCelli <= newCelli)
434  {
435  nbrCelli = -1;
436  }
437  else
438  {
439  // A region boundary? - treat like an external face
440  label ownRegion = cellToRegion[oldCelli];
441  label neiRegion = cellToRegion[cellOrder[nbrCelli]];
442 
443  if (ownRegion != neiRegion)
444  {
445  nbrCelli = -1;
446  }
447  }
448  }
449 
450  nbr.push_back(nbrCelli);
451  }
452 
453  Foam::sortedOrder(nbr, order);
454 
455  for (const label index : order)
456  {
457  if (nbr[index] >= 0)
458  {
459  oldToNewFace[cFaces[index]] = newFacei++;
460  }
461  }
462  }
463 
464  // This seems to be broken...
465 
466  // Do region interfaces
467  {
468  const label nRegions = Foam::max(cellToRegion)+1;
469 
470  // Sort in increasing region
472 
473  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
474  {
475  label ownRegion = cellToRegion[mesh.faceOwner()[facei]];
476  label neiRegion = cellToRegion[mesh.faceNeighbour()[facei]];
477 
478  if (ownRegion != neiRegion)
479  {
480  sortKey[facei] =
481  (
482  Foam::min(ownRegion, neiRegion)*nRegions
483  + Foam::max(ownRegion, neiRegion)
484  );
485  }
486  }
487 
488  sortKey.sort();
489 
490  // Extract.
491  forAll(sortKey, i)
492  {
493  label key = sortKey[i];
494 
495  if (key == labelMax)
496  {
497  break;
498  }
499 
500  oldToNewFace[sortKey.indices()[i]] = newFacei++;
501  }
502  }
503 
504  // Leave patch faces intact.
505  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
506  {
507  oldToNewFace[facei] = facei;
508  }
509 
510 
511  // Check done all faces.
512  forAll(oldToNewFace, facei)
513  {
514  if (oldToNewFace[facei] == -1)
515  {
517  << "Did not determine new position for face " << facei
518  << abort(FatalError);
519  }
520  }
521 
522  return invert(mesh.nFaces(), oldToNewFace);
523 }
524 
525 
526 // cellOrder: old cell for every new cell
527 // faceOrder: old face for every new face. Ordering of boundary faces not
528 // changed.
529 autoPtr<mapPolyMesh> reorderMesh
530 (
531  polyMesh& mesh,
532  const labelList& cellOrder,
533  const labelList& faceOrder
534 )
535 {
536  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
537  labelList reverseFaceOrder(invert(faceOrder.size(), faceOrder));
538 
539  faceList newFaces(reorder(reverseFaceOrder, mesh.faces()));
540  labelList newOwner
541  (
542  renumber
543  (
544  reverseCellOrder,
545  reorder(reverseFaceOrder, mesh.faceOwner())
546  )
547  );
548  labelList newNeighbour
549  (
550  renumber
551  (
552  reverseCellOrder,
553  reorder(reverseFaceOrder, mesh.faceNeighbour())
554  )
555  );
556 
557  // Check if any faces need swapping.
558  labelHashSet flipFaceFlux(newOwner.size());
559  forAll(newNeighbour, facei)
560  {
561  if (newOwner[facei] > newNeighbour[facei])
562  {
563  std::swap(newOwner[facei], newNeighbour[facei]);
564  newFaces[facei].flip();
565  flipFaceFlux.insert(facei);
566  }
567  }
568 
570  labelList patchSizes(patches.size());
571  labelList patchStarts(patches.size());
572  labelList oldPatchNMeshPoints(patches.size());
573  labelListList patchPointMap(patches.size());
574 
575  forAll(patches, patchi)
576  {
577  patchSizes[patchi] = patches[patchi].size();
578  patchStarts[patchi] = patches[patchi].start();
579  oldPatchNMeshPoints[patchi] = patches[patchi].nPoints();
580  patchPointMap[patchi] = identity(patches[patchi].nPoints());
581  }
582 
584  (
585  autoPtr<pointField>(), // <- null: leaves points untouched
586  autoPtr<faceList>::New(std::move(newFaces)),
587  autoPtr<labelList>::New(std::move(newOwner)),
588  autoPtr<labelList>::New(std::move(newNeighbour)),
589  patchSizes,
590  patchStarts,
591  true
592  );
593 
594 
595  // Re-do the faceZones
596  {
597  faceZoneMesh& faceZones = mesh.faceZones();
598  faceZones.clearAddressing();
599  forAll(faceZones, zoneI)
600  {
601  faceZone& fZone = faceZones[zoneI];
602  labelList newAddressing(fZone.size());
603  boolList newFlipMap(fZone.size());
604  forAll(fZone, i)
605  {
606  label oldFacei = fZone[i];
607  newAddressing[i] = reverseFaceOrder[oldFacei];
608  if (flipFaceFlux.found(newAddressing[i]))
609  {
610  newFlipMap[i] = !fZone.flipMap()[i];
611  }
612  else
613  {
614  newFlipMap[i] = fZone.flipMap()[i];
615  }
616  }
617  labelList newToOld(sortedOrder(newAddressing));
618  fZone.resetAddressing
619  (
620  labelUIndList(newAddressing, newToOld)(),
621  boolUIndList(newFlipMap, newToOld)()
622  );
623  }
624  }
625  // Re-do the cellZones
626  {
627  cellZoneMesh& cellZones = mesh.cellZones();
628  cellZones.clearAddressing();
629  forAll(cellZones, zoneI)
630  {
631  cellZones[zoneI] = labelUIndList
632  (
633  reverseCellOrder,
634  cellZones[zoneI]
635  )();
636  Foam::sort(cellZones[zoneI]);
637  }
638  }
639 
640 
642  (
643  mesh, // const polyMesh& mesh,
644  mesh.nPoints(), // nOldPoints,
645  mesh.nFaces(), // nOldFaces,
646  mesh.nCells(), // nOldCells,
647  identity(mesh.nPoints()), // pointMap,
648  List<objectMap>(), // pointsFromPoints,
649  faceOrder, // faceMap,
650  List<objectMap>(), // facesFromPoints,
651  List<objectMap>(), // facesFromEdges,
652  List<objectMap>(), // facesFromFaces,
653  cellOrder, // cellMap,
654  List<objectMap>(), // cellsFromPoints,
655  List<objectMap>(), // cellsFromEdges,
656  List<objectMap>(), // cellsFromFaces,
657  List<objectMap>(), // cellsFromCells,
658  identity(mesh.nPoints()), // reversePointMap,
659  reverseFaceOrder, // reverseFaceMap,
660  reverseCellOrder, // reverseCellMap,
661  flipFaceFlux, // flipFaceFlux,
662  patchPointMap, // patchPointMap,
663  labelListList(), // pointZoneMap,
664  labelListList(), // faceZonePointMap,
665  labelListList(), // faceZoneFaceMap,
666  labelListList(), // cellZoneMap,
667  pointField(), // preMotionPoints,
668  patchStarts, // oldPatchStarts,
669  oldPatchNMeshPoints, // oldPatchNMeshPoints
670  autoPtr<scalarField>() // oldCellVolumes
671  );
672 }
673 
674 
675 // Return new to old cell numbering, region-wise
676 CompactListList<label> regionRenumber
677 (
678  const renumberMethod& method,
679  const fvMesh& mesh,
680  const labelList& cellToRegion,
681  label nRegions = -1 // Number of regions or auto-detect
682 )
683 {
684  // Info<< "Determining cell order:" << endl;
685 
686  if (nRegions < 0)
687  {
688  nRegions = Foam::max(cellToRegion)+1;
689  }
690 
691  // Initially the per-region cell selection
692  CompactListList<label> regionCellOrder
693  (
694  invertOneToManyCompact(nRegions, cellToRegion)
695  );
696 
697  if (method.no_topology())
698  {
699  // Special case when renumberMesh is only used for decomposition.
700  // - can skip generating the connectivity
701  // - nonetheless calculate the order in case it is non-identity
702 
703  timer.resetTimeIncrement();
704 
705  forAll(regionCellOrder, regioni)
706  {
707  // Note: cellMap is identical to regionToCells[regioni]
708  // since it is already sorted
709 
710  labelList subCellOrder =
711  method.renumber(regionCellOrder[regioni].size());
712 
713  // Per region reordering (inplace but with SubList)
714  regionCellOrder[regioni] =
715  labelUIndList(regionCellOrder[regioni], subCellOrder)();
716  }
717 
718  timings[TimingType::RENUMBER] += timer.timeIncrement();
719  }
720  else if (method.needs_mesh())
721  {
722  timer.resetTimeIncrement();
723 
724  forAll(regionCellOrder, regioni)
725  {
726  // Info<< " region " << regioni << " starts at "
727  // << regionCellOrder.localStart(regioni) << nl;
728 
729  // No parallel communication
730  const bool oldParRun = UPstream::parRun(false);
731 
732  // Get local sub-mesh
733  fvMeshSubset subsetter(mesh, regionCellOrder[regioni]);
734 
735  // Note: cellMap will be identical to regionToCells[regioni]
736  // (assuming they are properly sorted!)
737  const labelList& cellMap = subsetter.cellMap();
738 
739  labelList subCellOrder = method.renumber(subsetter.subMesh());
740 
741  UPstream::parRun(oldParRun); // Restore parallel state
742 
743  // Per region reordering
744  regionCellOrder[regioni] = labelUIndList(cellMap, subCellOrder);
745  }
746 
747  timings[TimingType::RENUMBER] += timer.timeIncrement();
748  }
749  else
750  {
751  timer.resetTimeIncrement();
752 
753  // Create adjacency matrix of the full mesh and subset subsequently.
754  // This is more efficient than creating adjacency matrices of
755  // sub-meshes.
756 
757  // No parallel communication
758  const bool oldParRun = UPstream::parRun(false);
759 
760  // The local connectivity of the full (non-subsetted) mesh
761  CompactListList<label> meshCellCells;
762  globalMeshData::calcCellCells(mesh, meshCellCells);
763  UPstream::parRun(oldParRun); // Restore parallel state
764 
765  timings[TimingType::CELL_CELLS] += timer.timeIncrement();
766 
767  // For the respective subMesh selections
768  bitSet subsetCells(mesh.nCells());
769 
770  forAll(regionCellOrder, regioni)
771  {
772  // Info<< " region " << regioni << " starts at "
773  // << regionCellOrder.localStart(regioni) << nl;
774 
775  subsetCells = false;
776  subsetCells.set(regionCellOrder[regioni]);
777 
778  // Connectivity of local sub-mesh
779  labelList cellMap;
780  CompactListList<label> subCellCells =
781  subsetAdjacency(subsetCells, meshCellCells, cellMap);
782 
783  timings[TimingType::CELL_CELLS] += timer.timeIncrement();
784 
785  // No parallel communication
786  const bool oldParRun = UPstream::parRun(false);
787 
788  labelList subCellOrder = method.renumber(subCellCells);
789 
790  UPstream::parRun(oldParRun); // Restore parallel state
791 
792  // Per region reordering
793  regionCellOrder[regioni] = labelUIndList(cellMap, subCellOrder);
794 
795  timings[TimingType::RENUMBER] += timer.timeIncrement();
796  }
797  }
798  // Info<< endl;
799 
800  return regionCellOrder;
801 }
802 
803 
804 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
805 
806 int main(int argc, char *argv[])
807 {
808  argList::noFunctionObjects(); // Never use function objects
809 
810  timeSelector::addOptions_singleTime(); // Single-time options
811 
813  (
814  "Renumber mesh cells to reduce the bandwidth. Use the -lib option or"
815  " dictionary 'libs' entry to load additional libraries"
816  );
817 
819  (
820  "Test without writing. "
821  "Changes -write-maps to write VTK output."
822  );
824 
825  argList::addOption("dict", "file", "Alternative renumberMeshDict");
826 
828  (
829  "frontWidth",
830  "Calculate the RMS of the front-width"
831  );
832 
834  (
835  "decompose",
836  "Aggregate initially with a decomposition method (serial only)"
837  );
838 
840  (
841  "write-maps",
842  "Write renumber mappings"
843  );
844 
846  (
847  "no-fields",
848  "Suppress renumbering of fields (eg, when they are only uniform)"
849  );
850 
852  (
853  "list-renumber",
854  "List available renumbering methods"
855  );
856 
858  (
859  "renumber-method",
860  "name",
861  "Specify renumber method (default: CuthillMcKee) without dictionary"
862  );
863 
865  (
866  "renumber-coeffs",
867  "string-content",
868  "Specify renumber coefficients (dictionary content) as string. "
869  "eg, 'reverse true;'"
870  );
871 
872  #include "addAllRegionOptions.H"
873  #include "addOverwriteOption.H"
874 
875  // -------------------------
876 
877  #include "setRootCase.H"
878 
879  {
880  bool listOptions = false;
881 
882  if (args.found("list-renumber"))
883  {
884  listOptions = true;
885  Info<< nl
886  << "Available renumber methods:" << nl
887  << " "
889  << nl;
890  }
891 
892  if (listOptions)
893  {
894  return 0;
895  }
896  }
897 
898  const bool dryrun = args.dryRun();
899 
900  const bool readDict = args.found("dict");
901  const bool doDecompose = args.found("decompose");
902  const bool overwrite = args.found("overwrite");
903  const bool doFields = !args.found("no-fields");
904  const bool doFrontWidth = args.found("frontWidth") && !doDecompose;
905 
906  word renumberMethodName;
907  args.readIfPresent("renumber-method", renumberMethodName);
908 
909  if (doDecompose && UPstream::parRun())
910  {
912  << "Cannot use -decompose option in parallel ... giving up" << nl
913  << exit(FatalError);
914  }
915 
916  #include "createTime.H"
917 
918  // Allow override of decomposeParDict location
919  fileName decompDictFile;
920  if
921  (
922  args.readIfPresent("decomposeParDict", decompDictFile)
923  && !decompDictFile.empty() && !decompDictFile.isAbsolute()
924  )
925  {
926  decompDictFile = runTime.globalPath()/decompDictFile;
927  }
928 
929  // Get region names
930  #include "getAllRegionOptions.H"
931 
932  // Set time from specified time options, or force start from Time=0
934 
935  // Capture current time information for non-overwrite
937  (
940  );
941 
942  // Start/reset all timings
943  timer.resetTime();
944  timings = Foam::zero{};
945 
946  #include "createNamedMeshes.H"
947 
948  timings[TimingType::READ_MESH] += timer.timeIncrement();
949 
950 
951  for (fvMesh& mesh : meshes)
952  {
953  // Reset time in case of multiple meshes and not overwrite
954  if (!overwrite)
955  {
956  runTime.setTime(startTime.first(), startTime.second());
957  }
958 
959  const word oldInstance = mesh.pointsInstance();
960 
961  label band;
962  scalar profile;
963  scalar sumSqrIntersect;
964  getBand
965  (
966  doFrontWidth,
967  mesh.nCells(),
968  mesh.faceOwner(),
970  band,
971  profile,
972  sumSqrIntersect
973  );
974 
975  reduce(band, maxOp<label>());
976  reduce(profile, sumOp<scalar>());
977 
978  Info<< "Mesh " << mesh.name()
979  << " size: " << mesh.globalData().nTotalCells() << nl
980  << "Before renumbering" << nl
981  << " band : " << band << nl
982  << " profile : " << profile << nl;
983 
984  if (doFrontWidth)
985  {
986  reduce(sumSqrIntersect, sumOp<scalar>());
987  scalar rmsFrontwidth = Foam::sqrt
988  (
989  sumSqrIntersect/mesh.globalData().nTotalCells()
990  );
991 
992  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
993  }
994 
995  Info<< endl;
996 
997  bool sortCoupledFaceCells = false;
998  bool writeMaps = args.found("write-maps");
999  bool orderPoints = false;
1000  bool useRegionFaceOrder = false;
1001  label blockSize = 0;
1002 
1003  // Construct renumberMethod
1004  dictionary renumberDict;
1005  autoPtr<renumberMethod> renumberPtr;
1006 
1007  if (readDict)
1008  {
1009  const word dictName("renumberMeshDict");
1010  #include "setSystemMeshDictionaryIO.H"
1011 
1012  Info<< "Renumber according to " << dictIO.name() << nl << endl;
1013 
1014  renumberDict = IOdictionary::readContents(dictIO);
1015 
1016  renumberPtr = renumberMethod::New(renumberDict);
1017 
1018  // This options are likely orthogonal to -decompose mode
1019  if (!doDecompose)
1020  {
1021  sortCoupledFaceCells =
1022  renumberDict.getOrDefault("sortCoupledFaceCells", false);
1023 
1024  if (sortCoupledFaceCells)
1025  {
1026  Info<< "Sorting cells on coupled boundaries to be last."
1027  << nl << endl;
1028  }
1029 
1030  blockSize = renumberDict.getOrDefault<label>("blockSize", 0);
1031 
1032  if (blockSize > 0)
1033  {
1034  Info<< "Ordering cells into regions of size " << blockSize
1035  << " (using decomposition);"
1036  << " ordering faces into region-internal"
1037  << " and region-external."
1038  << nl << endl;
1039  }
1040 
1041  if (blockSize > 0)
1042  {
1043  useRegionFaceOrder =
1044  renumberDict.getOrDefault("regionFaceOrder", false);
1045  }
1046  }
1047 
1048  orderPoints = renumberDict.getOrDefault("orderPoints", false);
1049  if (orderPoints)
1050  {
1051  Info<< "Ordering points into internal and boundary points."
1052  << nl << endl;
1053  }
1054 
1055  if
1056  (
1057  renumberDict.readIfPresent("writeMaps", writeMaps)
1058  && writeMaps
1059  )
1060  {
1061  Info<< "Writing renumber maps (new to old) to polyMesh."
1062  << nl << endl;
1063  }
1064  }
1065  else
1066  {
1067  if (args.found("renumber-coeffs"))
1068  {
1069  ITstream is = args.lookup("renumber-coeffs");
1070 
1071  is >> renumberDict;
1072 
1073  Info<< "Specified renumber coefficients:" << nl
1074  << renumberDict << nl;
1075  }
1076 
1077  if (!renumberMethodName.empty())
1078  {
1079  // Specify/override the "method" within dictionary
1080  renumberDict.set("method", renumberMethodName);
1081  renumberPtr = renumberMethod::New(renumberDict);
1082  }
1083  else if (renumberDict.found("method"))
1084  {
1085  // Use the "method" type within dictionary
1086  renumberPtr = renumberMethod::New(renumberDict);
1087  }
1088  }
1089 
1090  // Default renumbering method
1091  if (!renumberPtr)
1092  {
1093  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
1094  Info<< "Using renumber-method: " << renumberPtr().type()
1095  << " [default]" << endl;
1096  }
1097  else
1098  {
1099  Info<< "Using renumber-method: " << renumberPtr().type()
1100  << endl;
1101  }
1102 
1103 
1104  // Read parallel reconstruct maps
1105  labelIOList cellProcAddressing
1106  (
1107  IOobject
1108  (
1109  "cellProcAddressing",
1110  mesh.facesInstance(),
1112  mesh,
1115  ),
1116  labelList()
1117  );
1118 
1119  labelIOList faceProcAddressing
1120  (
1121  IOobject
1122  (
1123  "faceProcAddressing",
1124  mesh.facesInstance(),
1126  mesh,
1129  ),
1130  labelList()
1131  );
1132  labelIOList pointProcAddressing
1133  (
1134  IOobject
1135  (
1136  "pointProcAddressing",
1137  mesh.pointsInstance(),
1139  mesh,
1142  ),
1143  labelList()
1144  );
1145  labelIOList boundaryProcAddressing
1146  (
1147  IOobject
1148  (
1149  "boundaryProcAddressing",
1150  mesh.pointsInstance(),
1152  mesh,
1155  ),
1156  labelList()
1157  );
1158 
1159 
1160  // List of stored objects to clear from mesh (after reading)
1161  DynamicList<regIOobject*> storedObjects;
1162 
1163  if (!dryrun && doFields)
1164  {
1165  Info<< nl << "Reading fields" << nl;
1166 
1167  timer.resetTimeIncrement();
1168 
1169  IOobjectList objects(mesh, runTime.timeName());
1170  storedObjects.reserve(objects.size());
1171 
1172  const predicates::always nameMatcher;
1173 
1174  // Read GeometricFields
1175 
1176  #undef doLocalCode
1177  #define doLocalCode(FieldType) \
1178  readFields<FieldType>(mesh, objects, nameMatcher, storedObjects);
1179 
1180  // Read volume fields
1186 
1187  // Read internal fields
1193 
1194  // Read surface fields
1200 
1201 
1202  // Read point fields
1203  const pointMesh& pMesh = pointMesh::New(mesh);
1204  #undef doLocalCode
1205  #define doLocalCode(FieldType) \
1206  readFields<FieldType>(pMesh, objects, nameMatcher, storedObjects);
1207 
1213 
1214  #undef doLocalCode
1215 
1216  timings[TimingType::READ_FIELDS] += timer.timeIncrement();
1217 
1218  // Write loaded fields when mesh.write() is called
1219  for (auto* fldptr : storedObjects)
1220  {
1221  fldptr->writeOpt(IOobject::AUTO_WRITE);
1222  }
1223  }
1224 
1225 
1226  // Read sets
1227  PtrList<cellSet> cellSets;
1228  PtrList<faceSet> faceSets;
1229  PtrList<pointSet> pointSets;
1230 
1231  if (!dryrun)
1232  {
1233  // Read sets
1234  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1235  ReadFields(objects, cellSets);
1236  ReadFields(objects, faceSets);
1237  ReadFields(objects, pointSets);
1238  }
1239 
1240 
1241  Info<< endl;
1242 
1243  // From renumbering:
1244  // - from new cell/face back to original cell/face
1245  labelList cellOrder;
1246  labelList faceOrder;
1247 
1248  if (blockSize > 0 && !doDecompose)
1249  {
1250  timer.resetTimeIncrement();
1251 
1252  // Renumbering in two phases. Should be done in one so mapping of
1253  // fields is done correctly!
1254 
1255  label nBlocks = mesh.nCells()/blockSize;
1256  Info<< "nBlocks = " << nBlocks << endl;
1257 
1258  // Read decompositionMethod dictionary
1259  dictionary decomposeDict(renumberDict.subDict("blockCoeffs"));
1260  decomposeDict.set("numberOfSubdomains", nBlocks);
1261 
1262  // No parallel communication
1263  const bool oldParRun = UPstream::parRun(false);
1264 
1265  autoPtr<decompositionMethod> decomposePtr =
1266  decompositionMethod::New(decomposeDict);
1267 
1268  labelList cellToRegion
1269  (
1270  decomposePtr().decompose
1271  (
1272  mesh,
1273  mesh.cellCentres()
1274  )
1275  );
1276 
1277  UPstream::parRun(oldParRun); // Restore parallel state
1278 
1279  timings[TimingType::DECOMPOSE] += timer.timeIncrement();
1280 
1281  // For debugging: write out region
1282  createScalarField
1283  (
1284  mesh,
1285  "cellDist",
1286  cellToRegion
1287  )().write();
1288 
1289  Info<< nl << "Written decomposition as volScalarField to "
1290  << "cellDist for use in postprocessing."
1291  << nl << endl;
1292 
1293 
1294  CompactListList<label> regionCellOrder =
1295  regionRenumber
1296  (
1297  renumberPtr(),
1298  mesh,
1299  cellToRegion,
1300  decomposePtr().nDomains()
1301  );
1302 
1303  cellOrder = regionCellOrder.values();
1304 
1305  // Determine new to old face order with new cell numbering
1306  if (useRegionFaceOrder)
1307  {
1308  faceOrder = getRegionFaceOrder(mesh, cellOrder, cellToRegion);
1309  }
1310  else
1311  {
1312  faceOrder = getFaceOrder(mesh, cellOrder);
1313  }
1314  }
1315  else
1316  {
1317  if (doDecompose)
1318  {
1319  // Two-step renumbering.
1320  // 1. decompose into regions (like decomposePar)
1321  // 2. renumber each sub-region
1322 
1323  timer.resetTimeIncrement();
1324 
1325  // Read decompositionMethod dictionary
1326  IOdictionary decomposeDict
1327  (
1329  (
1330  IOobject
1331  (
1333  runTime.time().system(),
1334  mesh.regionName(), // region (if non-default)
1335  runTime,
1339  ),
1340  decompDictFile
1341  )
1342  );
1343 
1344  // No parallel communication
1345  const bool oldParRun = UPstream::parRun(false);
1346 
1347  autoPtr<decompositionMethod> decomposePtr
1348  = decompositionMethod::New(decomposeDict);
1349 
1350  labelList cellToRegion
1351  (
1352  decomposePtr().decompose
1353  (
1354  mesh,
1355  mesh.cellCentres()
1356  )
1357  );
1358 
1359  timings[TimingType::DECOMPOSE] += timer.timeIncrement();
1360 
1361  UPstream::parRun(oldParRun); // Restore parallel state
1362 
1363  CompactListList<label> regionCellOrder =
1364  regionRenumber
1365  (
1366  renumberPtr(),
1367  mesh,
1368  cellToRegion,
1369  decomposePtr().nDomains()
1370  );
1371 
1372  cellOrder = regionCellOrder.values();
1373 
1374  // HACK - retain partitioning information for possible use
1375  // at a later stage
1376  mesh.data().set
1377  (
1378  "requested.partition-offsets",
1379  regionCellOrder.offsets()
1380  );
1381  }
1382  else
1383  {
1384  // Determines sorted back to original cell ordering
1385 
1386  const auto& method = renumberPtr();
1387 
1388  timer.resetTimeIncrement();
1389 
1390  if (method.no_topology())
1391  {
1392  cellOrder = method.renumber(mesh.nCells());
1393  }
1394  else
1395  {
1396  cellOrder = method.renumber(mesh);
1397  }
1398 
1399  timings[TimingType::RENUMBER] += timer.timeIncrement();
1400  }
1401 
1402 
1403  if (sortCoupledFaceCells)
1404  {
1405  // Change order so all coupled patch faceCells are at the end.
1407 
1408  // Collect all boundary cells on coupled patches
1409  label nBndCells = 0;
1410  forAll(pbm, patchi)
1411  {
1412  if (pbm[patchi].coupled())
1413  {
1414  nBndCells += pbm[patchi].size();
1415  }
1416  }
1417 
1418  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
1419 
1420  labelList bndCells(nBndCells);
1421  labelList bndCellMap(nBndCells);
1422  nBndCells = 0;
1423  forAll(pbm, patchi)
1424  {
1425  if (pbm[patchi].coupled())
1426  {
1427  const labelUList& faceCells = pbm[patchi].faceCells();
1428  forAll(faceCells, i)
1429  {
1430  label celli = faceCells[i];
1431 
1432  if (reverseCellOrder[celli] != -1)
1433  {
1434  bndCells[nBndCells] = celli;
1435  bndCellMap[nBndCells++] =
1436  reverseCellOrder[celli];
1437  reverseCellOrder[celli] = -1;
1438  }
1439  }
1440  }
1441  }
1442  bndCells.resize(nBndCells);
1443  bndCellMap.resize(nBndCells);
1444 
1445  // Sort
1446  labelList order(Foam::sortedOrder(bndCellMap));
1447 
1448  // Redo newReverseCellOrder
1449  labelList newReverseCellOrder(mesh.nCells(), -1);
1450 
1451  label sortedI = mesh.nCells();
1452  forAllReverse(order, i)
1453  {
1454  label origCelli = bndCells[order[i]];
1455  newReverseCellOrder[origCelli] = --sortedI;
1456  }
1457 
1458  Info<< "Ordered all " << nBndCells
1459  << " cells with a coupled face"
1460  << " to the end of the cell list, starting at " << sortedI
1461  << endl;
1462 
1463  // Compact
1464  sortedI = 0;
1465  forAll(cellOrder, newCelli)
1466  {
1467  label origCelli = cellOrder[newCelli];
1468  if (newReverseCellOrder[origCelli] == -1)
1469  {
1470  newReverseCellOrder[origCelli] = sortedI++;
1471  }
1472  }
1473 
1474  // Update sorted back to original (unsorted) map
1475  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
1476  }
1477 
1478 
1479  // Determine new to old face order with new cell numbering
1480  faceOrder = getFaceOrder(mesh, cellOrder);
1481  }
1482 
1483 
1484  if (!overwrite)
1485  {
1486  ++runTime;
1487  }
1488 
1489 
1490  // Change the mesh.
1491  autoPtr<mapPolyMesh> map = reorderMesh(mesh, cellOrder, faceOrder);
1492 
1493 
1494  if (orderPoints)
1495  {
1496  polyTopoChange meshMod(mesh);
1497  autoPtr<mapPolyMesh> pointOrderMap = meshMod.changeMesh
1498  (
1499  mesh,
1500  false, // inflate
1501  true, // syncParallel
1502  false, // orderCells
1503  orderPoints // orderPoints
1504  );
1505 
1506  // Combine point reordering into map.
1507  const_cast<labelList&>(map().pointMap()) = labelUIndList
1508  (
1509  map().pointMap(),
1510  pointOrderMap().pointMap()
1511  )();
1512 
1514  (
1515  pointOrderMap().reversePointMap(),
1516  const_cast<labelList&>(map().reversePointMap())
1517  );
1518  }
1519 
1520 
1521  // Update fields
1522  mesh.updateMesh(map());
1523 
1524  // Update proc maps
1525  if (cellProcAddressing.headerOk())
1526  {
1527  if (returnReduceAnd(cellProcAddressing.size() == mesh.nCells()))
1528  {
1529  Info<< "Renumbering processor cell decomposition map "
1530  << cellProcAddressing.name() << endl;
1531 
1532  cellProcAddressing = labelList
1533  (
1534  labelUIndList(cellProcAddressing, map().cellMap())
1535  );
1536  }
1537  else
1538  {
1539  Info<< "Not writing inconsistent processor cell decomposition"
1540  << " map " << cellProcAddressing.filePath() << endl;
1541  cellProcAddressing.writeOpt(IOobject::NO_WRITE);
1542  }
1543  }
1544  else
1545  {
1546  cellProcAddressing.writeOpt(IOobject::NO_WRITE);
1547  }
1548 
1549  if (faceProcAddressing.headerOk())
1550  {
1551  if (returnReduceAnd(faceProcAddressing.size() == mesh.nFaces()))
1552  {
1553  Info<< "Renumbering processor face decomposition map "
1554  << faceProcAddressing.name() << endl;
1555 
1556  faceProcAddressing = labelList
1557  (
1558  labelUIndList(faceProcAddressing, map().faceMap())
1559  );
1560 
1561  // Detect any flips.
1562  const labelHashSet& fff = map().flipFaceFlux();
1563  for (const label facei : fff)
1564  {
1565  label masterFacei = faceProcAddressing[facei];
1566 
1567  faceProcAddressing[facei] = -masterFacei;
1568 
1569  if (masterFacei == 0)
1570  {
1572  << " masterFacei:" << masterFacei
1573  << exit(FatalError);
1574  }
1575  }
1576  }
1577  else
1578  {
1579  Info<< "Not writing inconsistent processor face decomposition"
1580  << " map " << faceProcAddressing.filePath() << endl;
1581  faceProcAddressing.writeOpt(IOobject::NO_WRITE);
1582  }
1583  }
1584  else
1585  {
1586  faceProcAddressing.writeOpt(IOobject::NO_WRITE);
1587  }
1588 
1589  if (pointProcAddressing.headerOk())
1590  {
1591  if (returnReduceAnd(pointProcAddressing.size() == mesh.nPoints()))
1592  {
1593  Info<< "Renumbering processor point decomposition map "
1594  << pointProcAddressing.name() << endl;
1595 
1596  pointProcAddressing = labelList
1597  (
1598  labelUIndList(pointProcAddressing, map().pointMap())
1599  );
1600  }
1601  else
1602  {
1603  Info<< "Not writing inconsistent processor point decomposition"
1604  << " map " << pointProcAddressing.filePath() << endl;
1605  pointProcAddressing.writeOpt(IOobject::NO_WRITE);
1606  }
1607  }
1608  else
1609  {
1610  pointProcAddressing.writeOpt(IOobject::NO_WRITE);
1611  }
1612 
1613  if (boundaryProcAddressing.headerOk())
1614  {
1615  if
1616  (
1618  (
1619  boundaryProcAddressing.size() == mesh.boundaryMesh().size()
1620  )
1621  )
1622  {
1623  // No renumbering needed
1624  }
1625  else
1626  {
1627  Info<< "Not writing inconsistent processor patch decomposition"
1628  << " map " << boundaryProcAddressing.filePath() << endl;
1629  boundaryProcAddressing.writeOpt(IOobject::NO_WRITE);
1630  }
1631  }
1632  else
1633  {
1634  boundaryProcAddressing.writeOpt(IOobject::NO_WRITE);
1635  }
1636 
1637 
1638 
1639 
1640  // Move mesh (since morphing might not do this)
1641  if (map().hasMotionPoints())
1642  {
1643  mesh.movePoints(map().preMotionPoints());
1644  }
1645 
1646 
1647  {
1648  label band;
1649  scalar profile;
1650  scalar sumSqrIntersect;
1651  getBand
1652  (
1653  doFrontWidth,
1654  mesh.nCells(),
1655  mesh.faceOwner(),
1656  mesh.faceNeighbour(),
1657  band,
1658  profile,
1659  sumSqrIntersect
1660  );
1661  reduce(band, maxOp<label>());
1662  reduce(profile, sumOp<scalar>());
1663 
1664  Info<< "After renumbering";
1665  if (doDecompose)
1666  {
1667  Info<< " [values are misleading with -decompose option]";
1668  }
1669 
1670  Info<< nl
1671  << " band : " << band << nl
1672  << " profile : " << profile << nl;
1673 
1674  if (doFrontWidth)
1675  {
1676  reduce(sumSqrIntersect, sumOp<scalar>());
1677  scalar rmsFrontwidth = Foam::sqrt
1678  (
1679  sumSqrIntersect/mesh.globalData().nTotalCells()
1680  );
1681 
1682  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1683  }
1684 
1685  Info<< endl;
1686  }
1687 
1688  if (orderPoints)
1689  {
1690  // Force edge calculation (since only reason that points would
1691  // need to be sorted)
1692  (void)mesh.edges();
1693 
1694  label nTotPoints = returnReduce
1695  (
1696  mesh.nPoints(),
1697  sumOp<label>()
1698  );
1699  label nTotIntPoints = returnReduce
1700  (
1702  sumOp<label>()
1703  );
1704 
1705  label nTotEdges = returnReduce
1706  (
1707  mesh.nEdges(),
1708  sumOp<label>()
1709  );
1710  label nTotIntEdges = returnReduce
1711  (
1712  mesh.nInternalEdges(),
1713  sumOp<label>()
1714  );
1715  label nTotInt0Edges = returnReduce
1716  (
1718  sumOp<label>()
1719  );
1720  label nTotInt1Edges = returnReduce
1721  (
1723  sumOp<label>()
1724  );
1725 
1726  Info<< "Points:" << nl
1727  << " total : " << nTotPoints << nl
1728  << " internal: " << nTotIntPoints << nl
1729  << " boundary: " << nTotPoints-nTotIntPoints << nl
1730  << "Edges:" << nl
1731  << " total : " << nTotEdges << nl
1732  << " internal: " << nTotIntEdges << nl
1733  << " internal using 0 boundary points: "
1734  << nTotInt0Edges << nl
1735  << " internal using 1 boundary points: "
1736  << nTotInt1Edges-nTotInt0Edges << nl
1737  << " internal using 2 boundary points: "
1738  << nTotIntEdges-nTotInt1Edges << nl
1739  << " boundary: " << nTotEdges-nTotIntEdges << nl
1740  << endl;
1741  }
1742 
1743 
1744  if (dryrun)
1745  {
1746  if (writeMaps)
1747  {
1748  fileName file = mesh.time().globalPath()/"renumberMesh";
1749 
1750  const word& regionDir = mesh.regionName();
1751 
1752  if (!regionDir.empty())
1753  {
1754  file += "_" + regionDir;
1755  }
1756 
1757  // VTK output
1758  {
1759  const vtk::vtuCells cells(mesh);
1760 
1762  (
1763  mesh,
1764  cells,
1765  file,
1767  );
1768 
1771  writer.writeCellIDs();
1772 
1773  labelList ids;
1774 
1775  if (UPstream::parRun())
1776  {
1777  const label cellOffset =
1779 
1780  ids.resize(mesh.nCells());
1782  (
1783  map().cellMap().cbegin(),
1784  map().cellMap().cend(),
1785  ids.begin(),
1786  [=](const label id) { return (id + cellOffset); }
1787  );
1788 
1789  writer.writeCellData("origCellID", ids);
1790  writer.writeProcIDs();
1791  }
1792  else
1793  {
1794  writer.writeCellData("origCellID", map().cellMap());
1795 
1796  // Recover any partition information (in serial)
1797  globalIndex partitionOffsets;
1798  if
1799  (
1801  (
1802  "requested.partition-offsets",
1803  partitionOffsets.offsets()
1804  )
1805  )
1806  {
1807  if (partitionOffsets.totalSize() != mesh.nCells())
1808  {
1810  << "Requested partition total-size: "
1811  << partitionOffsets.totalSize()
1812  << " mesh total-size: "
1813  << mesh.nCells()
1814  << " ... ignoring" << endl;
1815 
1816  partitionOffsets.clear();
1817  }
1818  }
1819 
1820  ids.resize(partitionOffsets.totalSize());
1821 
1822  for (const label proci : partitionOffsets.allProcs())
1823  {
1824  ids.slice(partitionOffsets.range(proci)) = proci;
1825  }
1826 
1827  if (!partitionOffsets.empty())
1828  {
1829  writer.writeCellData("procID", ids);
1830  }
1831  }
1832 
1833  Info<< "Wrote renumbered mesh to "
1835  << " for visualization."
1836  << endl;
1837  }
1838  }
1839  }
1840  else
1841  {
1842  timer.resetTimeIncrement();
1843 
1844  if (overwrite)
1845  {
1846  mesh.setInstance(oldInstance);
1847  }
1848  else
1849  {
1851  }
1852 
1853  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1854 
1855  // Remove old procAddressing files
1857 
1858  // Update refinement data
1859 
1860  hexRef8Data refData
1861  (
1862  IOobject
1863  (
1864  "dummy",
1865  mesh.facesInstance(),
1867  mesh,
1871  )
1872  );
1873  refData.updateMesh(map());
1874  refData.write();
1875 
1876  // Update sets
1877  topoSet::updateMesh(mesh.facesInstance(), map(), cellSets);
1878  topoSet::updateMesh(mesh.facesInstance(), map(), faceSets);
1879  topoSet::updateMesh(mesh.facesInstance(), map(), pointSets);
1880 
1881  mesh.write();
1882 
1883  timings[TimingType::WRITING] += timer.timeIncrement();
1884 
1885  if (writeMaps)
1886  {
1887  // For debugging: write out region
1888  createScalarField
1889  (
1890  mesh,
1891  "origCellID",
1892  map().cellMap()
1893  )().write();
1894 
1895  createScalarField
1896  (
1897  mesh,
1898  "cellID",
1899  identity(mesh.nCells())
1900  )().write();
1901 
1902  Info<< nl
1903  << "Wrote current cellID and origCellID as volScalarField"
1904  << " for use in postprocessing." << nl << endl;
1905 
1906  IOobject meshMapIO
1907  (
1908  "map-name",
1909  mesh.facesInstance(),
1911  mesh,
1915  );
1916 
1917  meshMapIO.resetHeader("cellMap");
1918  IOList<label>::writeContents(meshMapIO, map().cellMap());
1919 
1920  meshMapIO.resetHeader("faceMap");
1921  IOList<label>::writeContents(meshMapIO, map().faceMap());
1922 
1923  meshMapIO.resetHeader("pointMap");
1924  IOList<label>::writeContents(meshMapIO, map().pointMap());
1925  }
1926  }
1927 
1928  // Cleanup loaded fields
1929  while (!storedObjects.empty())
1930  {
1931  storedObjects.back()->checkOut();
1932  storedObjects.pop_back();
1933  }
1934  }
1935 
1936  Info<< nl
1937  << "Timings:" << nl
1938  << " read mesh : " << timings[TimingType::READ_MESH] << nl
1939  << " read fields : " << timings[TimingType::READ_FIELDS] << nl
1940  << " decompose : " << timings[TimingType::DECOMPOSE] << nl
1941  << " cell-cells : " << timings[TimingType::CELL_CELLS] << nl
1942  << " renumber : " << timings[TimingType::RENUMBER] << nl
1943  << " write : " << timings[TimingType::WRITING] << nl
1944  << "TotalTime = " << timer.elapsedTime() << " s" << nl
1945  << nl;
1946 
1948 
1949  Info<< "End\n" << endl;
1950 
1951  return 0;
1952 }
1953 
1954 
1955 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:202
Foam::surfaceFields.
const polyBoundaryMesh & pbm
CompactListList< label > subsetAdjacency(const bitSet &select, const CompactListList< label > &input, labelList &subMap)
const Type & value() const noexcept
Return const reference to value.
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:562
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:477
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
const labelList & offsets() const noexcept
Return the offset table (= size()+1)
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:944
A class for handling file names.
Definition: fileName.H:72
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
label size() const noexcept
The primary size (the number of rows/sublists)
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:844
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:153
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:130
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:600
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:56
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1107
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:205
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
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 within a list.
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Implements a timeout mechanism via sigalarm.
Definition: timer.H:82
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
virtual bool no_topology() const
Renumbering method without mesh or cell-cell topology (very special case)
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:411
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:519
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:884
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1681
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath...
Definition: TimePathsI.H:122
virtual const meshState & data() const noexcept
Const reference to the mesh and solver state data.
Definition: polyMesh.H:559
bool writeProcIDs()
Write processor ids for each line as CellData or for each point as PointData, depending on isPointDat...
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:389
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:689
static const word canonicalName
The canonical name ("decomposeParDict") under which the MeshObject is registered. ...
fileName globalPath() const
The global path for the case = rootPath/globalCaseName.
Definition: TimePathsI.H:108
static bool setTimeIfPresent(Time &runTime, const argList &args, const bool forceInitial=false)
Set the runTime based on -constant (if present), -time (value), or -latestTime.
Definition: timeSelector.C:314
const cellList & cells() const
Field reading functions for post-processing utilities.
static bool isAbsolute(const std::string &str)
Return true if filename starts with a &#39;/&#39; or &#39;\&#39; or (windows-only) with a filesystem-root.
Definition: fileNameI.H:129
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
static void writeContents(const IOobject &io, const UList< T > &content)
Write contents. The IOobject is never registered.
Definition: IOList.C:194
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const globalIndex & globalMeshCellAddr() const noexcept
Global numbering for mesh cells.
label nFaces() const noexcept
Number of mesh faces.
void push_back(const T &val)
Append an element at the end of the list.
Definition: ListI.H:221
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
#define doLocalCode(FieldType, Variable)
static dictionary readContents(const IOobject &io)
Read and return contents, testing for "dictionary" type. The IOobject will not be registered...
Definition: IOdictionary.C:89
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
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)
void reduce(T &value, [[maybe_unused]] BinaryOp bop, [[maybe_unused]] const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce)
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:299
const labelList & offsets() const noexcept
Const-access to the offsets.
Definition: globalIndexI.H:252
label nInternal1Edges() const
Internal edges using 0 or 1 boundary point.
void writeCellData(const word &fieldName, const UList< Type > &field)
Write primitive field of CellData.
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
virtual labelList renumber(const label nCells) const
Return the cell visit order (from ordered back to original cell id) based solely on the number of cel...
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
Abstract base class for renumbering.
Unary and binary predicates that always return true, useful for templating.
Definition: predicates.H:53
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition: globalIndex.H:76
CompactListList< label > invertOneToManyCompact(const label len, const labelUList &map)
Invert one-to-many compact map. Unmapped elements will be size 0.
Definition: ListOps.C:175
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
label totalSize() const noexcept
The total addressed size, which corresponds to the end offset and also the sum of all localSizes...
Definition: globalIndexI.H:218
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:838
labelRange range(label proci) const noexcept
Return start/size range of proci data.
Definition: globalIndexI.H:332
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:960
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
SubList< T > slice(const label pos, label len=-1)
Return SubList slice (non-const access) - no range checking.
Definition: SubList.H:258
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
virtual bool writeGeometry()
Write patch topology.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:133
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
static void addDryRunOption(const string &usage, bool advanced=false)
Enable a &#39;dry-run&#39; bool option, with usage information.
Definition: argList.C:519
label nPoints
const Time & time() const noexcept
Return time registry.
const word & executable() const noexcept
Name of executable without the path.
Definition: argListI.H:44
virtual bool beginCellData(label nFields=0)
Begin CellData output section for specified number of fields.
int dryRun() const noexcept
Return the dry-run flag.
Definition: argListI.H:109
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:283
Cuthill-McKee renumbering (CM or RCM)
static void addOptions_singleTime()
Add single-time timeSelector options to argList::validOptions()
Definition: timeSelector.C:146
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
label nInternalPoints() const noexcept
Points not on boundary.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const word & system() const noexcept
Return system name.
Definition: TimePathsI.H:137
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1101
Reading is optional [identical to LAZY_READ].
static const word null
An empty word.
Definition: word.H:84
virtual void resetAddressing(faceZone &&zn)
Move reset addressing and flip map from another zone.
Definition: faceZone.C:500
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:400
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1296
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:905
label nInternalFaces() const noexcept
Number of internal faces.
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1088
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
void updateMesh(const mapPolyMesh &)
Update local numbering for changed mesh.
Definition: hexRef8Data.C:308
errorManip< error > abort(error &err)
Definition: errorManip.H:139
label nEdges() const
Number of mesh edges.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:410
ITstream lookup(const word &optName) const
Return an input stream from the named option.
Definition: argListI.H:177
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1068
virtual bool needs_mesh() const
Renumbering method requires a polyMesh for its topology.
label nInternalEdges() const
Internal edges using 0,1 or 2 boundary points.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh...
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...
A packed storage of objects of type <T> using an offset table for access.
static wordList supportedMethods()
Return a list of the known methods.
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:535
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
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
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:671
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:608
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
Total global number of mesh cells.
labelRange allProcs() const noexcept
Range of process indices for all addressed offsets (processes)
Definition: globalIndexI.H:236
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
Definition: foamGltfBase.H:105
Write an OpenFOAM volume (internal) geometry and internal fields as a vtu file or a legacy vtk file...
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:28
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name...
Definition: instant.H:53
static autoPtr< renumberMethod > New(const dictionary &dict)
Construct/select a renumbering method.
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 within a list.
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:58
#define WarningInFunction
Report a warning using Foam::Warning.
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.
const word & regionDir
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: PtrList.H:56
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:52
Direct mesh changes based on v1.3 polyTopoChange syntax.
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition: HashTable.C:729
const polyBoundaryMesh & patches
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Nothing to be read.
Required Classes.
Automatically write from objectRegistry::writeObject()
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:679
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:389
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
label localStart(const label proci) const
Start of proci data.
Definition: globalIndexI.H:283
bool empty() const noexcept
Check for default constructed or total-size == 0.
Definition: globalIndexI.H:173
const List< T > & values() const noexcept
Return the packed values.
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:75
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:765
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:315
static void calcCellCells(const polyMesh &mesh, const labelUList &agglom, const label nLocalCoarse, const bool parallel, CompactListList< label > &cellCells)
Determine (local or global) cellCells from mesh agglomeration.
bool readIfPresent(const word &optName, T &val) const
Read a value from the named option if present.
Definition: argListI.H:316
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
const fileName & output() const noexcept
The current output file name.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
bool coupled
Foam::argList args(argc, argv)
void resetHeader(const word &newName=word::null)
Clear various bits (headerClassName, note, sizeof...) that would be obtained when reading from a file...
Definition: IOobject.C:644
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
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:188
virtual void updateMesh(const mapPolyMesh &morphMap)
Update any stored data for new labels. Not implemented.
Definition: topoSet.C:687
Foam::label startTime
label nInternal0Edges() const
Internal edges (i.e. not on boundary face) using no boundary point.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:592
static void writeMaps(Ostream &os, const word &key, const labelListList &maps)
static const word & regionName(const word &region)
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:796
Do not request registration (bool: false)
A deep-copy description of an OpenFOAM volume mesh in data structures suitable for VTK UnstructuredGr...
Definition: foamVtuCells.H:73
void clear()
Reset to be empty (no offsets)
Definition: globalIndexI.H:264
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
An input stream of tokens.
Definition: ITstream.H:52
Namespace for OpenFOAM.
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...
static IOobject selectIO(const IOobject &io, const fileName &altFile, const word &ioName="")
Return the IOobject, but also consider an alternative file name.
Definition: IOobject.C:256
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:217
Starts timing and returns elapsed time from start. Uses std::chrono::high_resolution_clock for better...
Definition: clockTime.H:58
Required Variables.