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-2024 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 = 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  min(ownRegion, neiRegion)*nRegions
482  +max(ownRegion, neiRegion);
483  }
484  }
485 
486  sortKey.sort();
487 
488  // Extract.
489  forAll(sortKey, i)
490  {
491  label key = sortKey[i];
492 
493  if (key == labelMax)
494  {
495  break;
496  }
497 
498  oldToNewFace[sortKey.indices()[i]] = newFacei++;
499  }
500  }
501 
502  // Leave patch faces intact.
503  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
504  {
505  oldToNewFace[facei] = facei;
506  }
507 
508 
509  // Check done all faces.
510  forAll(oldToNewFace, facei)
511  {
512  if (oldToNewFace[facei] == -1)
513  {
515  << "Did not determine new position for face " << facei
516  << abort(FatalError);
517  }
518  }
519 
520  return invert(mesh.nFaces(), oldToNewFace);
521 }
522 
523 
524 // cellOrder: old cell for every new cell
525 // faceOrder: old face for every new face. Ordering of boundary faces not
526 // changed.
527 autoPtr<mapPolyMesh> reorderMesh
528 (
529  polyMesh& mesh,
530  const labelList& cellOrder,
531  const labelList& faceOrder
532 )
533 {
534  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
535  labelList reverseFaceOrder(invert(faceOrder.size(), faceOrder));
536 
537  faceList newFaces(reorder(reverseFaceOrder, mesh.faces()));
538  labelList newOwner
539  (
540  renumber
541  (
542  reverseCellOrder,
543  reorder(reverseFaceOrder, mesh.faceOwner())
544  )
545  );
546  labelList newNeighbour
547  (
548  renumber
549  (
550  reverseCellOrder,
551  reorder(reverseFaceOrder, mesh.faceNeighbour())
552  )
553  );
554 
555  // Check if any faces need swapping.
556  labelHashSet flipFaceFlux(newOwner.size());
557  forAll(newNeighbour, facei)
558  {
559  if (newOwner[facei] > newNeighbour[facei])
560  {
561  std::swap(newOwner[facei], newNeighbour[facei]);
562  newFaces[facei].flip();
563  flipFaceFlux.insert(facei);
564  }
565  }
566 
568  labelList patchSizes(patches.size());
569  labelList patchStarts(patches.size());
570  labelList oldPatchNMeshPoints(patches.size());
571  labelListList patchPointMap(patches.size());
572 
573  forAll(patches, patchi)
574  {
575  patchSizes[patchi] = patches[patchi].size();
576  patchStarts[patchi] = patches[patchi].start();
577  oldPatchNMeshPoints[patchi] = patches[patchi].nPoints();
578  patchPointMap[patchi] = identity(patches[patchi].nPoints());
579  }
580 
582  (
583  autoPtr<pointField>(), // <- null: leaves points untouched
584  autoPtr<faceList>::New(std::move(newFaces)),
585  autoPtr<labelList>::New(std::move(newOwner)),
586  autoPtr<labelList>::New(std::move(newNeighbour)),
587  patchSizes,
588  patchStarts,
589  true
590  );
591 
592 
593  // Re-do the faceZones
594  {
595  faceZoneMesh& faceZones = mesh.faceZones();
596  faceZones.clearAddressing();
597  forAll(faceZones, zoneI)
598  {
599  faceZone& fZone = faceZones[zoneI];
600  labelList newAddressing(fZone.size());
601  boolList newFlipMap(fZone.size());
602  forAll(fZone, i)
603  {
604  label oldFacei = fZone[i];
605  newAddressing[i] = reverseFaceOrder[oldFacei];
606  if (flipFaceFlux.found(newAddressing[i]))
607  {
608  newFlipMap[i] = !fZone.flipMap()[i];
609  }
610  else
611  {
612  newFlipMap[i] = fZone.flipMap()[i];
613  }
614  }
615  labelList newToOld(sortedOrder(newAddressing));
616  fZone.resetAddressing
617  (
618  labelUIndList(newAddressing, newToOld)(),
619  boolUIndList(newFlipMap, newToOld)()
620  );
621  }
622  }
623  // Re-do the cellZones
624  {
625  cellZoneMesh& cellZones = mesh.cellZones();
626  cellZones.clearAddressing();
627  forAll(cellZones, zoneI)
628  {
629  cellZones[zoneI] = labelUIndList
630  (
631  reverseCellOrder,
632  cellZones[zoneI]
633  )();
634  Foam::sort(cellZones[zoneI]);
635  }
636  }
637 
638 
640  (
641  mesh, // const polyMesh& mesh,
642  mesh.nPoints(), // nOldPoints,
643  mesh.nFaces(), // nOldFaces,
644  mesh.nCells(), // nOldCells,
645  identity(mesh.nPoints()), // pointMap,
646  List<objectMap>(), // pointsFromPoints,
647  faceOrder, // faceMap,
648  List<objectMap>(), // facesFromPoints,
649  List<objectMap>(), // facesFromEdges,
650  List<objectMap>(), // facesFromFaces,
651  cellOrder, // cellMap,
652  List<objectMap>(), // cellsFromPoints,
653  List<objectMap>(), // cellsFromEdges,
654  List<objectMap>(), // cellsFromFaces,
655  List<objectMap>(), // cellsFromCells,
656  identity(mesh.nPoints()), // reversePointMap,
657  reverseFaceOrder, // reverseFaceMap,
658  reverseCellOrder, // reverseCellMap,
659  flipFaceFlux, // flipFaceFlux,
660  patchPointMap, // patchPointMap,
661  labelListList(), // pointZoneMap,
662  labelListList(), // faceZonePointMap,
663  labelListList(), // faceZoneFaceMap,
664  labelListList(), // cellZoneMap,
665  pointField(), // preMotionPoints,
666  patchStarts, // oldPatchStarts,
667  oldPatchNMeshPoints, // oldPatchNMeshPoints
668  autoPtr<scalarField>() // oldCellVolumes
669  );
670 }
671 
672 
673 // Return new to old cell numbering, region-wise
674 CompactListList<label> regionRenumber
675 (
676  const renumberMethod& method,
677  const fvMesh& mesh,
678  const labelList& cellToRegion,
679  label nRegions = -1 // Number of regions or auto-detect
680 )
681 {
682  // Info<< "Determining cell order:" << endl;
683 
684  if (nRegions < 0)
685  {
686  nRegions = Foam::max(cellToRegion)+1;
687  }
688 
689  // Initially the per-region cell selection
690  CompactListList<label> regionCellOrder
691  (
692  invertOneToManyCompact(nRegions, cellToRegion)
693  );
694 
695  if (method.no_topology())
696  {
697  // Special case when renumberMesh is only used for decomposition.
698  // - can skip generating the connectivity
699  // - nonetheless calculate the order in case it is non-identity
700 
701  timer.resetTimeIncrement();
702 
703  forAll(regionCellOrder, regioni)
704  {
705  // Note: cellMap is identical to regionToCells[regioni]
706  // since it is already sorted
707 
708  labelList subCellOrder =
709  method.renumber(regionCellOrder[regioni].size());
710 
711  // Per region reordering (inplace but with SubList)
712  regionCellOrder[regioni] =
713  labelUIndList(regionCellOrder[regioni], subCellOrder)();
714  }
715 
716  timings[TimingType::RENUMBER] += timer.timeIncrement();
717  }
718  else if (method.needs_mesh())
719  {
720  timer.resetTimeIncrement();
721 
722  forAll(regionCellOrder, regioni)
723  {
724  // Info<< " region " << regioni << " starts at "
725  // << regionCellOrder.localStart(regioni) << nl;
726 
727  // No parallel communication
728  const bool oldParRun = UPstream::parRun(false);
729 
730  // Get local sub-mesh
731  fvMeshSubset subsetter(mesh, regionCellOrder[regioni]);
732 
733  // Note: cellMap will be identical to regionToCells[regioni]
734  // (assuming they are properly sorted!)
735  const labelList& cellMap = subsetter.cellMap();
736 
737  labelList subCellOrder = method.renumber(subsetter.subMesh());
738 
739  UPstream::parRun(oldParRun); // Restore parallel state
740 
741  // Per region reordering
742  regionCellOrder[regioni] = labelUIndList(cellMap, subCellOrder);
743  }
744 
745  timings[TimingType::RENUMBER] += timer.timeIncrement();
746  }
747  else
748  {
749  timer.resetTimeIncrement();
750 
751  // Create adjacency matrix of the full mesh and subset subsequently.
752  // This is more efficient than creating adjacency matrices of
753  // sub-meshes.
754 
755  // No parallel communication
756  const bool oldParRun = UPstream::parRun(false);
757 
758  // The local connectivity of the full (non-subsetted) mesh
759  CompactListList<label> meshCellCells;
760  globalMeshData::calcCellCells(mesh, meshCellCells);
761  UPstream::parRun(oldParRun); // Restore parallel state
762 
763  timings[TimingType::CELL_CELLS] += timer.timeIncrement();
764 
765  // For the respective subMesh selections
766  bitSet subsetCells(mesh.nCells());
767 
768  forAll(regionCellOrder, regioni)
769  {
770  // Info<< " region " << regioni << " starts at "
771  // << regionCellOrder.localStart(regioni) << nl;
772 
773  subsetCells = false;
774  subsetCells.set(regionCellOrder[regioni]);
775 
776  // Connectivity of local sub-mesh
777  labelList cellMap;
778  CompactListList<label> subCellCells =
779  subsetAdjacency(subsetCells, meshCellCells, cellMap);
780 
781  timings[TimingType::CELL_CELLS] += timer.timeIncrement();
782 
783  // No parallel communication
784  const bool oldParRun = UPstream::parRun(false);
785 
786  labelList subCellOrder = method.renumber(subCellCells);
787 
788  UPstream::parRun(oldParRun); // Restore parallel state
789 
790  // Per region reordering
791  regionCellOrder[regioni] = labelUIndList(cellMap, subCellOrder);
792 
793  timings[TimingType::RENUMBER] += timer.timeIncrement();
794  }
795  }
796  // Info<< endl;
797 
798  return regionCellOrder;
799 }
800 
801 
802 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
803 
804 int main(int argc, char *argv[])
805 {
806  argList::noFunctionObjects(); // Never use function objects
807 
808  timeSelector::addOptions_singleTime(); // Single-time options
809 
811  (
812  "Renumber mesh cells to reduce the bandwidth. Use the -lib option or"
813  " dictionary 'libs' entry to load additional libraries"
814  );
815 
817  (
818  "Test without writing. "
819  "Changes -write-maps to write VTK output."
820  );
822 
823  argList::addOption("dict", "file", "Alternative renumberMeshDict");
824 
826  (
827  "frontWidth",
828  "Calculate the RMS of the front-width"
829  );
830 
832  (
833  "decompose",
834  "Aggregate initially with a decomposition method (serial only)"
835  );
836 
838  (
839  "write-maps",
840  "Write renumber mappings"
841  );
842 
844  (
845  "no-fields",
846  "Suppress renumbering of fields (eg, when they are only uniform)"
847  );
848 
850  (
851  "list-renumber",
852  "List available renumbering methods"
853  );
854 
856  (
857  "renumber-method",
858  "name",
859  "Specify renumber method (default: CuthillMcKee) without dictionary"
860  );
861 
863  (
864  "renumber-coeffs",
865  "string-content",
866  "Specify renumber coefficients (dictionary content) as string. "
867  "eg, 'reverse true;'"
868  );
869 
870  #include "addAllRegionOptions.H"
871  #include "addOverwriteOption.H"
872 
873  // -------------------------
874 
875  #include "setRootCase.H"
876 
877  {
878  bool listOptions = false;
879 
880  if (args.found("list-renumber"))
881  {
882  listOptions = true;
883  Info<< nl
884  << "Available renumber methods:" << nl
885  << " "
887  << nl;
888  }
889 
890  if (listOptions)
891  {
892  return 0;
893  }
894  }
895 
896  const bool dryrun = args.dryRun();
897 
898  const bool readDict = args.found("dict");
899  const bool doDecompose = args.found("decompose");
900  const bool overwrite = args.found("overwrite");
901  const bool doFields = !args.found("no-fields");
902  const bool doFrontWidth = args.found("frontWidth") && !doDecompose;
903 
904  word renumberMethodName;
905  args.readIfPresent("renumber-method", renumberMethodName);
906 
907  if (doDecompose && UPstream::parRun())
908  {
910  << "Cannot use -decompose option in parallel ... giving up" << nl
911  << exit(FatalError);
912  }
913 
914  #include "createTime.H"
915 
916  // Allow override of decomposeParDict location
917  fileName decompDictFile;
918  if
919  (
920  args.readIfPresent("decomposeParDict", decompDictFile)
921  && !decompDictFile.empty() && !decompDictFile.isAbsolute()
922  )
923  {
924  decompDictFile = runTime.globalPath()/decompDictFile;
925  }
926 
927  // Get region names
928  #include "getAllRegionOptions.H"
929 
930  // Set time from specified time options, or force start from Time=0
932 
933  // Capture current time information for non-overwrite
935  (
938  );
939 
940  // Start/reset all timings
941  timer.resetTime();
942  timings = Foam::zero{};
943 
944  #include "createNamedMeshes.H"
945 
946  timings[TimingType::READ_MESH] += timer.timeIncrement();
947 
948 
949  for (fvMesh& mesh : meshes)
950  {
951  // Reset time in case of multiple meshes and not overwrite
952  if (!overwrite)
953  {
954  runTime.setTime(startTime.first(), startTime.second());
955  }
956 
957  const word oldInstance = mesh.pointsInstance();
958 
959  label band;
960  scalar profile;
961  scalar sumSqrIntersect;
962  getBand
963  (
964  doFrontWidth,
965  mesh.nCells(),
966  mesh.faceOwner(),
968  band,
969  profile,
970  sumSqrIntersect
971  );
972 
973  reduce(band, maxOp<label>());
974  reduce(profile, sumOp<scalar>());
975 
976  Info<< "Mesh " << mesh.name()
977  << " size: " << mesh.globalData().nTotalCells() << nl
978  << "Before renumbering" << nl
979  << " band : " << band << nl
980  << " profile : " << profile << nl;
981 
982  if (doFrontWidth)
983  {
984  reduce(sumSqrIntersect, sumOp<scalar>());
985  scalar rmsFrontwidth = Foam::sqrt
986  (
987  sumSqrIntersect/mesh.globalData().nTotalCells()
988  );
989 
990  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
991  }
992 
993  Info<< endl;
994 
995  bool sortCoupledFaceCells = false;
996  bool writeMaps = args.found("write-maps");
997  bool orderPoints = false;
998  bool useRegionFaceOrder = false;
999  label blockSize = 0;
1000 
1001  // Construct renumberMethod
1002  dictionary renumberDict;
1003  autoPtr<renumberMethod> renumberPtr;
1004 
1005  if (readDict)
1006  {
1007  const word dictName("renumberMeshDict");
1008  #include "setSystemMeshDictionaryIO.H"
1009 
1010  Info<< "Renumber according to " << dictIO.name() << nl << endl;
1011 
1012  renumberDict = IOdictionary::readContents(dictIO);
1013 
1014  renumberPtr = renumberMethod::New(renumberDict);
1015 
1016  // This options are likely orthogonal to -decompose mode
1017  if (!doDecompose)
1018  {
1019  sortCoupledFaceCells =
1020  renumberDict.getOrDefault("sortCoupledFaceCells", false);
1021 
1022  if (sortCoupledFaceCells)
1023  {
1024  Info<< "Sorting cells on coupled boundaries to be last."
1025  << nl << endl;
1026  }
1027 
1028  blockSize = renumberDict.getOrDefault<label>("blockSize", 0);
1029 
1030  if (blockSize > 0)
1031  {
1032  Info<< "Ordering cells into regions of size " << blockSize
1033  << " (using decomposition);"
1034  << " ordering faces into region-internal"
1035  << " and region-external."
1036  << nl << endl;
1037  }
1038 
1039  if (blockSize > 0)
1040  {
1041  useRegionFaceOrder =
1042  renumberDict.getOrDefault("regionFaceOrder", false);
1043  }
1044  }
1045 
1046  orderPoints = renumberDict.getOrDefault("orderPoints", false);
1047  if (orderPoints)
1048  {
1049  Info<< "Ordering points into internal and boundary points."
1050  << nl << endl;
1051  }
1052 
1053  if
1054  (
1055  renumberDict.readIfPresent("writeMaps", writeMaps)
1056  && writeMaps
1057  )
1058  {
1059  Info<< "Writing renumber maps (new to old) to polyMesh."
1060  << nl << endl;
1061  }
1062  }
1063  else
1064  {
1065  if (args.found("renumber-coeffs"))
1066  {
1067  ITstream is = args.lookup("renumber-coeffs");
1068 
1069  is >> renumberDict;
1070 
1071  Info<< "Specified renumber coefficients:" << nl
1072  << renumberDict << nl;
1073  }
1074 
1075  if (!renumberMethodName.empty())
1076  {
1077  // Specify/override the "method" within dictionary
1078  renumberDict.set("method", renumberMethodName);
1079  renumberPtr = renumberMethod::New(renumberDict);
1080  }
1081  else if (renumberDict.found("method"))
1082  {
1083  // Use the "method" type within dictionary
1084  renumberPtr = renumberMethod::New(renumberDict);
1085  }
1086  }
1087 
1088  // Default renumbering method
1089  if (!renumberPtr)
1090  {
1091  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
1092  Info<< "Using renumber-method: " << renumberPtr().type()
1093  << " [default]" << endl;
1094  }
1095  else
1096  {
1097  Info<< "Using renumber-method: " << renumberPtr().type()
1098  << endl;
1099  }
1100 
1101 
1102  // Read parallel reconstruct maps
1103  labelIOList cellProcAddressing
1104  (
1105  IOobject
1106  (
1107  "cellProcAddressing",
1108  mesh.facesInstance(),
1110  mesh,
1113  ),
1114  labelList()
1115  );
1116 
1117  labelIOList faceProcAddressing
1118  (
1119  IOobject
1120  (
1121  "faceProcAddressing",
1122  mesh.facesInstance(),
1124  mesh,
1127  ),
1128  labelList()
1129  );
1130  labelIOList pointProcAddressing
1131  (
1132  IOobject
1133  (
1134  "pointProcAddressing",
1135  mesh.pointsInstance(),
1137  mesh,
1140  ),
1141  labelList()
1142  );
1143  labelIOList boundaryProcAddressing
1144  (
1145  IOobject
1146  (
1147  "boundaryProcAddressing",
1148  mesh.pointsInstance(),
1150  mesh,
1153  ),
1154  labelList()
1155  );
1156 
1157 
1158  // List of stored objects to clear from mesh (after reading)
1159  DynamicList<regIOobject*> storedObjects;
1160 
1161  if (!dryrun && doFields)
1162  {
1163  Info<< nl << "Reading fields" << nl;
1164 
1165  timer.resetTimeIncrement();
1166 
1167  IOobjectList objects(mesh, runTime.timeName());
1168  storedObjects.reserve(objects.size());
1169 
1170  const predicates::always nameMatcher;
1171 
1172  // Read GeometricFields
1173 
1174  #undef doLocalCode
1175  #define doLocalCode(FieldType) \
1176  readFields<FieldType>(mesh, objects, nameMatcher, storedObjects);
1177 
1178  // Read volume fields
1184 
1185  // Read internal fields
1191 
1192  // Read surface fields
1198 
1199 
1200  // Read point fields
1201  const pointMesh& pMesh = pointMesh::New(mesh);
1202  #undef doLocalCode
1203  #define doLocalCode(FieldType) \
1204  readFields<FieldType>(pMesh, objects, nameMatcher, storedObjects);
1205 
1211 
1212  #undef doLocalCode
1213 
1214  timings[TimingType::READ_FIELDS] += timer.timeIncrement();
1215 
1216  // Write loaded fields when mesh.write() is called
1217  for (auto* fldptr : storedObjects)
1218  {
1219  fldptr->writeOpt(IOobject::AUTO_WRITE);
1220  }
1221  }
1222 
1223 
1224  // Read sets
1225  PtrList<cellSet> cellSets;
1226  PtrList<faceSet> faceSets;
1227  PtrList<pointSet> pointSets;
1228 
1229  if (!dryrun)
1230  {
1231  // Read sets
1232  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1233  ReadFields(objects, cellSets);
1234  ReadFields(objects, faceSets);
1235  ReadFields(objects, pointSets);
1236  }
1237 
1238 
1239  Info<< endl;
1240 
1241  // From renumbering:
1242  // - from new cell/face back to original cell/face
1243  labelList cellOrder;
1244  labelList faceOrder;
1245 
1246  if (blockSize > 0 && !doDecompose)
1247  {
1248  timer.resetTimeIncrement();
1249 
1250  // Renumbering in two phases. Should be done in one so mapping of
1251  // fields is done correctly!
1252 
1253  label nBlocks = mesh.nCells()/blockSize;
1254  Info<< "nBlocks = " << nBlocks << endl;
1255 
1256  // Read decompositionMethod dictionary
1257  dictionary decomposeDict(renumberDict.subDict("blockCoeffs"));
1258  decomposeDict.set("numberOfSubdomains", nBlocks);
1259 
1260  // No parallel communication
1261  const bool oldParRun = UPstream::parRun(false);
1262 
1263  autoPtr<decompositionMethod> decomposePtr =
1264  decompositionMethod::New(decomposeDict);
1265 
1266  labelList cellToRegion
1267  (
1268  decomposePtr().decompose
1269  (
1270  mesh,
1271  mesh.cellCentres()
1272  )
1273  );
1274 
1275  UPstream::parRun(oldParRun); // Restore parallel state
1276 
1277  timings[TimingType::DECOMPOSE] += timer.timeIncrement();
1278 
1279  // For debugging: write out region
1280  createScalarField
1281  (
1282  mesh,
1283  "cellDist",
1284  cellToRegion
1285  )().write();
1286 
1287  Info<< nl << "Written decomposition as volScalarField to "
1288  << "cellDist for use in postprocessing."
1289  << nl << endl;
1290 
1291 
1292  CompactListList<label> regionCellOrder =
1293  regionRenumber
1294  (
1295  renumberPtr(),
1296  mesh,
1297  cellToRegion,
1298  decomposePtr().nDomains()
1299  );
1300 
1301  cellOrder = regionCellOrder.values();
1302 
1303  // Determine new to old face order with new cell numbering
1304  if (useRegionFaceOrder)
1305  {
1306  faceOrder = getRegionFaceOrder(mesh, cellOrder, cellToRegion);
1307  }
1308  else
1309  {
1310  faceOrder = getFaceOrder(mesh, cellOrder);
1311  }
1312  }
1313  else
1314  {
1315  if (doDecompose)
1316  {
1317  // Two-step renumbering.
1318  // 1. decompose into regions (like decomposePar)
1319  // 2. renumber each sub-region
1320 
1321  timer.resetTimeIncrement();
1322 
1323  // Read decompositionMethod dictionary
1324  IOdictionary decomposeDict
1325  (
1327  (
1328  IOobject
1329  (
1331  runTime.time().system(),
1332  mesh.regionName(), // region (if non-default)
1333  runTime,
1337  ),
1338  decompDictFile
1339  )
1340  );
1341 
1342  // No parallel communication
1343  const bool oldParRun = UPstream::parRun(false);
1344 
1345  autoPtr<decompositionMethod> decomposePtr
1346  = decompositionMethod::New(decomposeDict);
1347 
1348  labelList cellToRegion
1349  (
1350  decomposePtr().decompose
1351  (
1352  mesh,
1353  mesh.cellCentres()
1354  )
1355  );
1356 
1357  timings[TimingType::DECOMPOSE] += timer.timeIncrement();
1358 
1359  UPstream::parRun(oldParRun); // Restore parallel state
1360 
1361  CompactListList<label> regionCellOrder =
1362  regionRenumber
1363  (
1364  renumberPtr(),
1365  mesh,
1366  cellToRegion,
1367  decomposePtr().nDomains()
1368  );
1369 
1370  cellOrder = regionCellOrder.values();
1371 
1372  // HACK - retain partitioning information for possible use
1373  // at a later stage
1374  mesh.data().set
1375  (
1376  "requested.partition-offsets",
1377  regionCellOrder.offsets()
1378  );
1379  }
1380  else
1381  {
1382  // Determines sorted back to original cell ordering
1383 
1384  const auto& method = renumberPtr();
1385 
1386  timer.resetTimeIncrement();
1387 
1388  if (method.no_topology())
1389  {
1390  cellOrder = method.renumber(mesh.nCells());
1391  }
1392  else
1393  {
1394  cellOrder = method.renumber(mesh);
1395  }
1396 
1397  timings[TimingType::RENUMBER] += timer.timeIncrement();
1398  }
1399 
1400 
1401  if (sortCoupledFaceCells)
1402  {
1403  // Change order so all coupled patch faceCells are at the end.
1405 
1406  // Collect all boundary cells on coupled patches
1407  label nBndCells = 0;
1408  forAll(pbm, patchi)
1409  {
1410  if (pbm[patchi].coupled())
1411  {
1412  nBndCells += pbm[patchi].size();
1413  }
1414  }
1415 
1416  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
1417 
1418  labelList bndCells(nBndCells);
1419  labelList bndCellMap(nBndCells);
1420  nBndCells = 0;
1421  forAll(pbm, patchi)
1422  {
1423  if (pbm[patchi].coupled())
1424  {
1425  const labelUList& faceCells = pbm[patchi].faceCells();
1426  forAll(faceCells, i)
1427  {
1428  label celli = faceCells[i];
1429 
1430  if (reverseCellOrder[celli] != -1)
1431  {
1432  bndCells[nBndCells] = celli;
1433  bndCellMap[nBndCells++] =
1434  reverseCellOrder[celli];
1435  reverseCellOrder[celli] = -1;
1436  }
1437  }
1438  }
1439  }
1440  bndCells.resize(nBndCells);
1441  bndCellMap.resize(nBndCells);
1442 
1443  // Sort
1444  labelList order(Foam::sortedOrder(bndCellMap));
1445 
1446  // Redo newReverseCellOrder
1447  labelList newReverseCellOrder(mesh.nCells(), -1);
1448 
1449  label sortedI = mesh.nCells();
1450  forAllReverse(order, i)
1451  {
1452  label origCelli = bndCells[order[i]];
1453  newReverseCellOrder[origCelli] = --sortedI;
1454  }
1455 
1456  Info<< "Ordered all " << nBndCells
1457  << " cells with a coupled face"
1458  << " to the end of the cell list, starting at " << sortedI
1459  << endl;
1460 
1461  // Compact
1462  sortedI = 0;
1463  forAll(cellOrder, newCelli)
1464  {
1465  label origCelli = cellOrder[newCelli];
1466  if (newReverseCellOrder[origCelli] == -1)
1467  {
1468  newReverseCellOrder[origCelli] = sortedI++;
1469  }
1470  }
1471 
1472  // Update sorted back to original (unsorted) map
1473  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
1474  }
1475 
1476 
1477  // Determine new to old face order with new cell numbering
1478  faceOrder = getFaceOrder(mesh, cellOrder);
1479  }
1480 
1481 
1482  if (!overwrite)
1483  {
1484  ++runTime;
1485  }
1486 
1487 
1488  // Change the mesh.
1489  autoPtr<mapPolyMesh> map = reorderMesh(mesh, cellOrder, faceOrder);
1490 
1491 
1492  if (orderPoints)
1493  {
1494  polyTopoChange meshMod(mesh);
1495  autoPtr<mapPolyMesh> pointOrderMap = meshMod.changeMesh
1496  (
1497  mesh,
1498  false, // inflate
1499  true, // syncParallel
1500  false, // orderCells
1501  orderPoints // orderPoints
1502  );
1503 
1504  // Combine point reordering into map.
1505  const_cast<labelList&>(map().pointMap()) = labelUIndList
1506  (
1507  map().pointMap(),
1508  pointOrderMap().pointMap()
1509  )();
1510 
1512  (
1513  pointOrderMap().reversePointMap(),
1514  const_cast<labelList&>(map().reversePointMap())
1515  );
1516  }
1517 
1518 
1519  // Update fields
1520  mesh.updateMesh(map());
1521 
1522  // Update proc maps
1523  if (cellProcAddressing.headerOk())
1524  {
1525  if (returnReduceAnd(cellProcAddressing.size() == mesh.nCells()))
1526  {
1527  Info<< "Renumbering processor cell decomposition map "
1528  << cellProcAddressing.name() << endl;
1529 
1530  cellProcAddressing = labelList
1531  (
1532  labelUIndList(cellProcAddressing, map().cellMap())
1533  );
1534  }
1535  else
1536  {
1537  Info<< "Not writing inconsistent processor cell decomposition"
1538  << " map " << cellProcAddressing.filePath() << endl;
1539  cellProcAddressing.writeOpt(IOobject::NO_WRITE);
1540  }
1541  }
1542  else
1543  {
1544  cellProcAddressing.writeOpt(IOobject::NO_WRITE);
1545  }
1546 
1547  if (faceProcAddressing.headerOk())
1548  {
1549  if (returnReduceAnd(faceProcAddressing.size() == mesh.nFaces()))
1550  {
1551  Info<< "Renumbering processor face decomposition map "
1552  << faceProcAddressing.name() << endl;
1553 
1554  faceProcAddressing = labelList
1555  (
1556  labelUIndList(faceProcAddressing, map().faceMap())
1557  );
1558 
1559  // Detect any flips.
1560  const labelHashSet& fff = map().flipFaceFlux();
1561  for (const label facei : fff)
1562  {
1563  label masterFacei = faceProcAddressing[facei];
1564 
1565  faceProcAddressing[facei] = -masterFacei;
1566 
1567  if (masterFacei == 0)
1568  {
1570  << " masterFacei:" << masterFacei
1571  << exit(FatalError);
1572  }
1573  }
1574  }
1575  else
1576  {
1577  Info<< "Not writing inconsistent processor face decomposition"
1578  << " map " << faceProcAddressing.filePath() << endl;
1579  faceProcAddressing.writeOpt(IOobject::NO_WRITE);
1580  }
1581  }
1582  else
1583  {
1584  faceProcAddressing.writeOpt(IOobject::NO_WRITE);
1585  }
1586 
1587  if (pointProcAddressing.headerOk())
1588  {
1589  if (returnReduceAnd(pointProcAddressing.size() == mesh.nPoints()))
1590  {
1591  Info<< "Renumbering processor point decomposition map "
1592  << pointProcAddressing.name() << endl;
1593 
1594  pointProcAddressing = labelList
1595  (
1596  labelUIndList(pointProcAddressing, map().pointMap())
1597  );
1598  }
1599  else
1600  {
1601  Info<< "Not writing inconsistent processor point decomposition"
1602  << " map " << pointProcAddressing.filePath() << endl;
1603  pointProcAddressing.writeOpt(IOobject::NO_WRITE);
1604  }
1605  }
1606  else
1607  {
1608  pointProcAddressing.writeOpt(IOobject::NO_WRITE);
1609  }
1610 
1611  if (boundaryProcAddressing.headerOk())
1612  {
1613  if
1614  (
1616  (
1617  boundaryProcAddressing.size() == mesh.boundaryMesh().size()
1618  )
1619  )
1620  {
1621  // No renumbering needed
1622  }
1623  else
1624  {
1625  Info<< "Not writing inconsistent processor patch decomposition"
1626  << " map " << boundaryProcAddressing.filePath() << endl;
1627  boundaryProcAddressing.writeOpt(IOobject::NO_WRITE);
1628  }
1629  }
1630  else
1631  {
1632  boundaryProcAddressing.writeOpt(IOobject::NO_WRITE);
1633  }
1634 
1635 
1636 
1637 
1638  // Move mesh (since morphing might not do this)
1639  if (map().hasMotionPoints())
1640  {
1641  mesh.movePoints(map().preMotionPoints());
1642  }
1643 
1644 
1645  {
1646  label band;
1647  scalar profile;
1648  scalar sumSqrIntersect;
1649  getBand
1650  (
1651  doFrontWidth,
1652  mesh.nCells(),
1653  mesh.faceOwner(),
1654  mesh.faceNeighbour(),
1655  band,
1656  profile,
1657  sumSqrIntersect
1658  );
1659  reduce(band, maxOp<label>());
1660  reduce(profile, sumOp<scalar>());
1661 
1662  Info<< "After renumbering";
1663  if (doDecompose)
1664  {
1665  Info<< " [values are misleading with -decompose option]";
1666  }
1667 
1668  Info<< nl
1669  << " band : " << band << nl
1670  << " profile : " << profile << nl;
1671 
1672  if (doFrontWidth)
1673  {
1674  reduce(sumSqrIntersect, sumOp<scalar>());
1675  scalar rmsFrontwidth = Foam::sqrt
1676  (
1677  sumSqrIntersect/mesh.globalData().nTotalCells()
1678  );
1679 
1680  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1681  }
1682 
1683  Info<< endl;
1684  }
1685 
1686  if (orderPoints)
1687  {
1688  // Force edge calculation (since only reason that points would
1689  // need to be sorted)
1690  (void)mesh.edges();
1691 
1692  label nTotPoints = returnReduce
1693  (
1694  mesh.nPoints(),
1695  sumOp<label>()
1696  );
1697  label nTotIntPoints = returnReduce
1698  (
1700  sumOp<label>()
1701  );
1702 
1703  label nTotEdges = returnReduce
1704  (
1705  mesh.nEdges(),
1706  sumOp<label>()
1707  );
1708  label nTotIntEdges = returnReduce
1709  (
1710  mesh.nInternalEdges(),
1711  sumOp<label>()
1712  );
1713  label nTotInt0Edges = returnReduce
1714  (
1716  sumOp<label>()
1717  );
1718  label nTotInt1Edges = returnReduce
1719  (
1721  sumOp<label>()
1722  );
1723 
1724  Info<< "Points:" << nl
1725  << " total : " << nTotPoints << nl
1726  << " internal: " << nTotIntPoints << nl
1727  << " boundary: " << nTotPoints-nTotIntPoints << nl
1728  << "Edges:" << nl
1729  << " total : " << nTotEdges << nl
1730  << " internal: " << nTotIntEdges << nl
1731  << " internal using 0 boundary points: "
1732  << nTotInt0Edges << nl
1733  << " internal using 1 boundary points: "
1734  << nTotInt1Edges-nTotInt0Edges << nl
1735  << " internal using 2 boundary points: "
1736  << nTotIntEdges-nTotInt1Edges << nl
1737  << " boundary: " << nTotEdges-nTotIntEdges << nl
1738  << endl;
1739  }
1740 
1741 
1742  if (dryrun)
1743  {
1744  if (writeMaps)
1745  {
1746  fileName file = mesh.time().globalPath()/"renumberMesh";
1747 
1748  const word& regionDir = mesh.regionName();
1749 
1750  if (!regionDir.empty())
1751  {
1752  file += "_" + regionDir;
1753  }
1754 
1755  // VTK output
1756  {
1757  const vtk::vtuCells cells(mesh);
1758 
1760  (
1761  mesh,
1762  cells,
1763  file,
1765  );
1766 
1769  writer.writeCellIDs();
1770 
1771  labelList ids;
1772 
1773  if (UPstream::parRun())
1774  {
1775  const label cellOffset =
1777 
1778  ids.resize(mesh.nCells());
1780  (
1781  map().cellMap().cbegin(),
1782  map().cellMap().cend(),
1783  ids.begin(),
1784  [=](const label id) { return (id + cellOffset); }
1785  );
1786 
1787  writer.writeCellData("origCellID", ids);
1788  writer.writeProcIDs();
1789  }
1790  else
1791  {
1792  writer.writeCellData("origCellID", map().cellMap());
1793 
1794  // Recover any partition information (in serial)
1795  globalIndex partitionOffsets;
1796  if
1797  (
1799  (
1800  "requested.partition-offsets",
1801  partitionOffsets.offsets()
1802  )
1803  )
1804  {
1805  if (partitionOffsets.totalSize() != mesh.nCells())
1806  {
1808  << "Requested partition total-size: "
1809  << partitionOffsets.totalSize()
1810  << " mesh total-size: "
1811  << mesh.nCells()
1812  << " ... ignoring" << endl;
1813 
1814  partitionOffsets.clear();
1815  }
1816  }
1817 
1818  ids.resize(partitionOffsets.totalSize());
1819 
1820  for (const label proci : partitionOffsets.allProcs())
1821  {
1822  ids.slice(partitionOffsets.range(proci)) = proci;
1823  }
1824 
1825  if (!partitionOffsets.empty())
1826  {
1827  writer.writeCellData("procID", ids);
1828  }
1829  }
1830 
1831  Info<< "Wrote renumbered mesh to "
1833  << " for visualization."
1834  << endl;
1835  }
1836  }
1837  }
1838  else
1839  {
1840  timer.resetTimeIncrement();
1841 
1842  if (overwrite)
1843  {
1844  mesh.setInstance(oldInstance);
1845  }
1846  else
1847  {
1849  }
1850 
1851  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1852 
1853  // Remove old procAddressing files
1855 
1856  // Update refinement data
1857 
1858  hexRef8Data refData
1859  (
1860  IOobject
1861  (
1862  "dummy",
1863  mesh.facesInstance(),
1865  mesh,
1869  )
1870  );
1871  refData.updateMesh(map());
1872  refData.write();
1873 
1874  // Update sets
1875  topoSet::updateMesh(mesh.facesInstance(), map(), cellSets);
1876  topoSet::updateMesh(mesh.facesInstance(), map(), faceSets);
1877  topoSet::updateMesh(mesh.facesInstance(), map(), pointSets);
1878 
1879  mesh.write();
1880 
1881  timings[TimingType::WRITING] += timer.timeIncrement();
1882 
1883  if (writeMaps)
1884  {
1885  // For debugging: write out region
1886  createScalarField
1887  (
1888  mesh,
1889  "origCellID",
1890  map().cellMap()
1891  )().write();
1892 
1893  createScalarField
1894  (
1895  mesh,
1896  "cellID",
1897  identity(mesh.nCells())
1898  )().write();
1899 
1900  Info<< nl
1901  << "Wrote current cellID and origCellID as volScalarField"
1902  << " for use in postprocessing." << nl << endl;
1903 
1904  IOobject meshMapIO
1905  (
1906  "map-name",
1907  mesh.facesInstance(),
1909  mesh,
1913  );
1914 
1915  meshMapIO.resetHeader("cellMap");
1916  IOListRef<label>(meshMapIO, map().cellMap()).write();
1917 
1918  meshMapIO.resetHeader("faceMap");
1919  IOListRef<label>(meshMapIO, map().faceMap()).write();
1920 
1921  meshMapIO.resetHeader("pointMap");
1922  IOListRef<label>(meshMapIO, map().pointMap()).write();
1923  }
1924  }
1925 
1926  // Cleanup loaded fields
1927  while (!storedObjects.empty())
1928  {
1929  storedObjects.back()->checkOut();
1930  storedObjects.pop_back();
1931  }
1932  }
1933 
1934  Info<< nl
1935  << "Timings:" << nl
1936  << " read mesh : " << timings[TimingType::READ_MESH] << nl
1937  << " read fields : " << timings[TimingType::READ_FIELDS] << nl
1938  << " decompose : " << timings[TimingType::DECOMPOSE] << nl
1939  << " cell-cells : " << timings[TimingType::CELL_CELLS] << nl
1940  << " renumber : " << timings[TimingType::RENUMBER] << nl
1941  << " write : " << timings[TimingType::WRITING] << nl
1942  << "TotalTime = " << timer.elapsedTime() << " s" << nl
1943  << nl;
1944 
1946 
1947  Info<< "End\n" << endl;
1948 
1949  return 0;
1950 }
1951 
1952 
1953 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
Foam::surfaceFields.
labelRange range(const label proci) const
Return start/size range of proci data.
Definition: globalIndexI.H:282
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: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
const labelList & offsets() const noexcept
Return the offset table (= size()+1)
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:905
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:859
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
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:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:54
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 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 within a list.
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:531
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:929
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
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:103
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: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
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
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:220
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
#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
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)
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:421
const labelList & offsets() const noexcept
Const-access to the offsets.
Definition: globalIndexI.H:205
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...
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.
Unary and binary predicates that always return true, useful for templating.
Definition: predicates.H:53
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:61
CompactListList< label > invertOneToManyCompact(const label len, const labelUList &map)
Invert one-to-many compact map. Unmapped elements will be size 0.
Definition: ListOps.C:176
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:165
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
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1005
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:250
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
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:130
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
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:504
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
int dryRun() const noexcept
Return the dry-run flag.
Definition: argListI.H:109
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
Cuthill-McKee renumbering (CM or RCM)
static void addOptions_singleTime()
Add single-time timeSelector options to argList::validOptions()
Definition: timeSelector.C:146
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:118
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].
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: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:902
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:1103
constexpr auto cend(const C &c) -> decltype(c.end())
Return const_iterator to the end of the container c.
Definition: stdFoam.H:223
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:392
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:1113
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 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...
A packed storage of objects of type <T> using an offset table for access.
virtual bool beginCellData(label nFields=0)
Begin CellData output section for specified number of fields.
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:520
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
static const word canonicalName
The canonical name ("decomposeParDict") under which the MeshObject is registered. ...
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:607
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:75
constexpr auto cbegin(const C &c) -> decltype(c.begin())
Return const_iterator to the beginning of the container c.
Definition: stdFoam.H:190
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:189
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))
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:30
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:59
#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
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: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:736
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:367
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool writeGeometry()
Write patch topology.
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:233
bool empty() const noexcept
Check for default constructed or total-size == 0.
Definition: globalIndexI.H:120
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:437
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
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::argList args(argc, argv)
A IOList wrapper for writing external data.
Definition: IOList.H:151
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:587
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:180
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:600
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:811
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:66
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
void clear()
Reset to be empty (no offsets)
Definition: globalIndexI.H:217
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.
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:225
Starts timing and returns elapsed time from start. Uses std::chrono::high_resolution_clock for better...
Definition: clockTime.H:58
Required Variables.