redistributePar.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-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2023 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  redistributePar
29 
30 Group
31  grpParallelUtilities
32 
33 Description
34  Redistributes existing decomposed mesh and fields according to the current
35  settings in the decomposeParDict file.
36 
37  Must be run on maximum number of source and destination processors.
38  Balances mesh and writes new mesh to new time directory.
39 
40  Can optionally run in decompose/reconstruct mode to decompose/reconstruct
41  mesh and fields.
42 
43 Usage
44  \b redistributePar [OPTION]
45 
46  Options:
47  - \par -decompose
48  Remove any existing \a processor subdirectories and decomposes the
49  mesh. Equivalent to running without processor subdirectories.
50 
51  - \par -reconstruct
52  Reconstruct mesh and fields (like reconstructParMesh+reconstructPar).
53 
54  - \par -newTimes
55  (in combination with -reconstruct) reconstruct only new times.
56 
57  - \par -dry-run
58  (not in combination with -reconstruct) Test without actually
59  decomposing.
60 
61  - \par -cellDist
62  not in combination with -reconstruct) Write the cell distribution
63  as a labelList, for use with 'manual'
64  decomposition method and as a volScalarField for visualization.
65 
66  - \par -region <regionName>
67  Distribute named region.
68 
69  - \par -allRegions
70  Distribute all regions in regionProperties. Does not check for
71  existence of processor*.
72 
73 \*---------------------------------------------------------------------------*/
74 
75 #include "argList.H"
76 #include "sigFpe.H"
77 #include "Time.H"
78 #include "fvMesh.H"
79 #include "fvMeshTools.H"
80 #include "fvMeshDistribute.H"
81 #include "fieldsDistributor.H"
82 #include "decompositionMethod.H"
83 #include "decompositionModel.H"
84 #include "timeSelector.H"
85 #include "PstreamReduceOps.H"
86 #include "volFields.H"
87 #include "surfaceFields.H"
89 #include "IOobjectList.H"
90 #include "globalIndex.H"
91 #include "loadOrCreateMesh.H"
92 #include "processorFvPatchField.H"
93 #include "topoSet.H"
94 #include "regionProperties.H"
95 
96 #include "parFvFieldDistributor.H"
98 #include "hexRef8Data.H"
99 #include "meshRefinement.H"
100 #include "pointFields.H"
101 
102 #include "faMeshSubset.H"
103 #include "faMeshTools.H"
104 #include "faMeshDistributor.H"
106 
107 #include "redistributeLagrangian.H"
108 
109 #include "cyclicACMIFvPatch.H"
111 #include "uncollatedFileOperation.H"
112 #include "collatedFileOperation.H"
113 
114 using namespace Foam;
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 const int debug(::Foam::debug::debugSwitch("redistributePar", 0));
119 
120 void createTimeDirs(const fileName& path)
121 {
122  // Get current set of local processor's time directories. Uses
123  // fileHandler
124  const instantList localTimeDirs(Time::findTimes(path, "constant"));
125 
126  instantList masterTimeDirs;
127  if (Pstream::master())
128  {
129  //const bool oldParRun = Pstream::parRun(false);
130  //timeDirs = Time::findTimes(path, "constant");
131  //Pstream::parRun(oldParRun); // Restore parallel state
132  masterTimeDirs = localTimeDirs;
133  }
134  Pstream::broadcast(masterTimeDirs);
135  //DebugVar(masterTimeDirs);
136  //DebugVar(localTimeDirs);
137 
138  // Sync any cached times (e.g. masterUncollatedFileOperation::times_)
139  // since only master would have done the findTimes
140  for (const instant& t : masterTimeDirs)
141  {
142  if (!localTimeDirs.found(t))
143  {
144  const fileName timePath(path/t.name());
145 
146  //Pout<< "Time:" << t << nl
147  // << " raw :" << timePath << nl
148  // << endl;
149  Foam::mkDir(timePath);
150  }
151  }
152 
153  // Just to make sure remove all state and re-scan
154  fileHandler().flush();
155  (void)Time::findTimes(path, "constant");
156 }
157 
158 
159 void copyUniform
160 (
161  const fileOperation& fh,
162  const bool decompose,
163  const bool reconstruct,
164  const word& readTimeName,
165  const objectRegistry& readDb,
166  const objectRegistry& writeDb
167 )
168 {
169  // Detect uniform/ at original database + time
170  const IOobject readIO("uniform", readTimeName, readDb);
171  const fileName readPath
172  (
173  fh.dirPath
174  (
175  false, // local directory
176  readIO,
177  false // do not search in time
178  )
179  );
180  //if (Pstream::master() && !readPath.empty())
181  if (!readPath.empty())
182  {
183  Info<< "Detected additional non-decomposed files in "
184  << readPath << endl;
185 
186  // readPath: searching is the same for all file handlers. Typical:
187  // <case>/0.1/uniform (parent dir, decompose mode)
188  // <case>/processor1/0.1/uniform (redistribute/reconstruct mode)
189  // <case>/processors2/0.1/uniform ,,
190  // writePath:
191  // uncollated : <case>/0.1/uniform (reconstruct mode). Should only
192  // be done by master
193  // uncollated : <case>/processorXXX/0.1/uniform. Should be done by all.
194  // collated : <case>/processors2/0.1/uniform. Should be done by
195  // local master only.
196 
197  // See what local directory
198  const IOobject writeIO("uniform", writeDb.time().timeName(), writeDb);
199  const fileName writePath
200  (
201  fh.objectPath
202  (
203  writeIO,
204  word::null
205  )
206  );
207  // Do we already have this directory?
208  const fileName currentPath(fh.dirPath(false, writeIO, false));
209 
210  if (::debug)
211  {
212  Pout<< " readPath :" << readPath << endl;
213  Pout<< " writePath :" << writePath << endl;
214  Pout<< " currentPath:" << currentPath << endl;
215  }
216 
217  if (readPath == writePath)
218  {
219  return;
220  }
221 
222  if (currentPath.empty())
223  {
224  if (decompose)
225  {
226  // All processors copy to destination
227  fh.cp(readPath, writePath);
228  }
229  else if (reconstruct)
230  {
231  // Only master
232  if (Pstream::master())
233  {
234  const bool oldParRun = Pstream::parRun(false);
235  fh.cp(readPath, writePath);
236  Pstream::parRun(oldParRun);
237  }
238  }
239  else
240  {
241  // Redistribute. If same destination path do only on master,
242  // if different path do on all processors. For now check
243  // if collated file handler only. tbd.
244  if (isA<fileOperations::collatedFileOperation>(fh))
245  {
246  // Collated
247  if (Pstream::master())
248  {
249  const bool oldParRun = Pstream::parRun(false);
250  fh.cp(readPath, writePath);
251  Pstream::parRun(oldParRun);
252  }
253  }
254  else
255  {
256  // Assume uncollated
257  fh.cp(readPath, writePath);
258  }
259  }
260  }
261  }
262 }
263 
264 
265 void printMeshData(const polyMesh& mesh)
266 {
267  // Collect all data on master
268 
269  labelListList patchNeiProcNo(Pstream::nProcs());
270  labelListList patchSize(Pstream::nProcs());
271  const labelList& pPatches = mesh.globalData().processorPatches();
272  patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
273  patchSize[Pstream::myProcNo()].setSize(pPatches.size());
274  forAll(pPatches, i)
275  {
276  const processorPolyPatch& ppp = refCast<const processorPolyPatch>
277  (
278  mesh.boundaryMesh()[pPatches[i]]
279  );
280  patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
281  patchSize[Pstream::myProcNo()][i] = ppp.size();
282  }
283  Pstream::gatherList(patchNeiProcNo);
284  Pstream::gatherList(patchSize);
285 
286 
287  // Print stats
288 
289  const globalIndex globalCells(mesh.nCells());
290  const globalIndex globalBoundaryFaces(mesh.nBoundaryFaces());
291 
292  label maxProcCells = 0;
293  label maxProcFaces = 0;
294  label totProcFaces = 0;
295  label maxProcPatches = 0;
296  label totProcPatches = 0;
297 
298  for (const int proci : Pstream::allProcs())
299  {
300  const label nLocalCells = globalCells.localSize(proci);
301  const label nBndFaces = globalBoundaryFaces.localSize(proci);
302 
303  Info<< nl
304  << "Processor " << proci;
305 
306  if (!nLocalCells)
307  {
308  Info<< " (empty)" << endl;
309  continue;
310  }
311  else
312  {
313  Info<< nl
314  << " Number of cells = " << nLocalCells << endl;
315  }
316 
317  label nProcFaces = 0;
318  const labelList& nei = patchNeiProcNo[proci];
319 
320  forAll(patchNeiProcNo[proci], i)
321  {
322  Info<< " Number of faces shared with processor "
323  << patchNeiProcNo[proci][i] << " = "
324  << patchSize[proci][i] << nl;
325 
326  nProcFaces += patchSize[proci][i];
327  }
328 
329  {
330  Info<< " Number of processor patches = " << nei.size() << nl
331  << " Number of processor faces = " << nProcFaces << nl
332  << " Number of boundary faces = "
333  << nBndFaces-nProcFaces << endl;
334  }
335 
336  maxProcCells = max(maxProcCells, nLocalCells);
337  totProcFaces += nProcFaces;
338  totProcPatches += nei.size();
339  maxProcFaces = max(maxProcFaces, nProcFaces);
340  maxProcPatches = max(maxProcPatches, nei.size());
341  }
342 
343  // Summary stats
344 
345  Info<< nl
346  << "Number of processor faces = " << (totProcFaces/2) << nl
347  << "Max number of cells = " << maxProcCells;
348 
349  if (maxProcCells != globalCells.totalSize())
350  {
351  scalar avgValue = scalar(globalCells.totalSize())/Pstream::nProcs();
352 
353  Info<< " (" << 100.0*(maxProcCells-avgValue)/avgValue
354  << "% above average " << avgValue << ')';
355  }
356  Info<< nl;
357 
358  Info<< "Max number of processor patches = " << maxProcPatches;
359  if (totProcPatches)
360  {
361  scalar avgValue = scalar(totProcPatches)/Pstream::nProcs();
362 
363  Info<< " (" << 100.0*(maxProcPatches-avgValue)/avgValue
364  << "% above average " << avgValue << ')';
365  }
366  Info<< nl;
367 
368  Info<< "Max number of faces between processors = " << maxProcFaces;
369  if (totProcFaces)
370  {
371  scalar avgValue = scalar(totProcFaces)/Pstream::nProcs();
372 
373  Info<< " (" << 100.0*(maxProcFaces-avgValue)/avgValue
374  << "% above average " << avgValue << ')';
375  }
376  Info<< nl << endl;
377 }
378 
379 
380 // Debugging: write volScalarField with decomposition for post processing.
381 void writeDecomposition
382 (
383  const word& name,
384  const fvMesh& mesh,
385  const labelUList& decomp
386 )
387 {
388  // Write the decomposition as labelList for use with 'manual'
389  // decomposition method.
390  labelIOList cellDecomposition
391  (
392  IOobject
393  (
394  "cellDecomposition",
395  mesh.facesInstance(), // mesh read from facesInstance
396  mesh,
400  ),
401  decomp
402  );
403  cellDecomposition.write();
404 
405  Info<< "Writing wanted cell distribution to volScalarField " << name
406  << " for postprocessing purposes." << nl << endl;
407 
408  volScalarField procCells
409  (
410  IOobject
411  (
412  name,
413  mesh.time().timeName(),
414  mesh,
418  ),
419  mesh,
422  );
423 
424  forAll(procCells, celli)
425  {
426  procCells[celli] = decomp[celli];
427  }
428 
429  procCells.correctBoundaryConditions();
430  procCells.write();
431 }
432 
433 
434 void determineDecomposition
435 (
436  const Time& baseRunTime,
437  const fileName& decompDictFile, // optional location for decomposeParDict
438  const bool decompose, // decompose, i.e. read from undecomposed case
439  const fileName& proc0CaseName,
440  const fvMesh& mesh,
441  const bool writeCellDist,
442 
443  label& nDestProcs,
444  labelList& decomp
445 )
446 {
447  // Read decomposeParDict (on all processors)
449  (
450  mesh,
451  decompDictFile
452  );
453 
454  decompositionMethod& decomposer = method.decomposer();
455 
456  if (!decomposer.parallelAware())
457  {
459  << "You have selected decomposition method \""
460  << decomposer.type() << "\n"
461  << " which does not synchronise decomposition across"
462  " processor patches.\n"
463  " You might want to select a decomposition method"
464  " that is aware of this. Continuing...." << endl;
465  }
466 
467  Time& tm = const_cast<Time&>(mesh.time());
468 
469  const bool oldProcCase = tm.processorCase();
470  if (Pstream::master() && decompose)
471  {
472  Info<< "Setting caseName to " << baseRunTime.caseName()
473  << " to read decomposeParDict" << endl;
474  tm.caseName() = baseRunTime.caseName();
475  tm.processorCase(false);
476  }
477 
478  scalarField cellWeights;
479  if (method.found("weightField"))
480  {
481  word weightName = method.get<word>("weightField");
482 
483  volScalarField weights
484  (
485  IOobject
486  (
487  weightName,
488  tm.timeName(),
489  mesh,
492  ),
493  mesh
494  );
495  cellWeights = weights.internalField();
496  }
497 
498  nDestProcs = decomposer.nDomains();
499  decomp = decomposer.decompose(mesh, cellWeights);
500 
501  if (Pstream::master() && decompose)
502  {
503  Info<< "Restoring caseName" << endl;
504  tm.caseName() = proc0CaseName;
505  tm.processorCase(oldProcCase);
506  }
507 
508  // Dump decomposition to volScalarField
509  if (writeCellDist)
510  {
511  // Note: on master make sure to write to processor0
512  if (decompose)
513  {
514  if (Pstream::master())
515  {
516  const bool oldParRun = Pstream::parRun(false);
517 
518  Info<< "Setting caseName to " << baseRunTime.caseName()
519  << " to write undecomposed cellDist" << endl;
520 
521  tm.caseName() = baseRunTime.caseName();
522  tm.processorCase(false);
523  writeDecomposition("cellDist", mesh, decomp);
524  Info<< "Restoring caseName" << endl;
525  tm.caseName() = proc0CaseName;
526  tm.processorCase(oldProcCase);
527 
528  Pstream::parRun(oldParRun);
529  }
530  }
531  else
532  {
533  writeDecomposition("cellDist", mesh, decomp);
534  }
535  }
536 }
537 
538 
539 // Variant of GeometricField::correctBoundaryConditions that only
540 // evaluates selected patch fields
541 template<class CoupledPatchType, class GeoField>
542 void correctCoupledBoundaryConditions(fvMesh& mesh)
543 {
544  for (GeoField& fld : mesh.sorted<GeoField>())
545  {
546  fld.boundaryFieldRef().template evaluateCoupled<CoupledPatchType>();
547  }
548 }
549 
550 
551 // Inplace redistribute mesh and any fields
552 autoPtr<mapDistributePolyMesh> redistributeAndWrite
553 (
554  refPtr<fileOperation>& writeHandler,
555  const Time& baseRunTime,
556  const fileName& proc0CaseName,
557 
558  // Controls
559  const bool doReadFields,
560  const bool decompose, // decompose, i.e. read from undecomposed case
561  const bool reconstruct,
562  const bool overwrite,
563 
564  // Decomposition information
565  const label nDestProcs,
566  const labelList& decomp,
567 
568  // Mesh information
569  const boolList& volMeshOnProc,
570  const fileName& volMeshInstance,
571  fvMesh& mesh
572 )
573 {
574  Time& runTime = const_cast<Time&>(mesh.time());
575  const bool oldProcCase = runTime.processorCase();
576 
578  //Info<< "Before distribution:" << endl;
579  //printMeshData(mesh);
580 
581  // Storage of fields
582 
583  PtrList<volScalarField> volScalarFields;
584  PtrList<volVectorField> volVectorFields;
585  PtrList<volSphericalTensorField> volSphereTensorFields;
586  PtrList<volSymmTensorField> volSymmTensorFields;
587  PtrList<volTensorField> volTensorFields;
588 
589  PtrList<surfaceScalarField> surfScalarFields;
590  PtrList<surfaceVectorField> surfVectorFields;
591  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
592  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
593  PtrList<surfaceTensorField> surfTensorFields;
594 
600 
601  PtrList<pointScalarField> pointScalarFields;
602  PtrList<pointVectorField> pointVectorFields;
603  PtrList<pointTensorField> pointTensorFields;
604  PtrList<pointSphericalTensorField> pointSphTensorFields;
605  PtrList<pointSymmTensorField> pointSymmTensorFields;
606 
607  // Self-contained pointMesh for reading pointFields
608  const pointMesh oldPointMesh(mesh);
609 
610  // Track how many (if any) pointFields are read/mapped
611  label nPointFields = 0;
612 
613  parPointFieldDistributor pointDistributor
614  (
615  oldPointMesh, // source mesh
616  false, // savePoints=false (ie, delay until later)
617  false // Do not write
618  );
619 
620 
621  if (doReadFields)
622  {
623  // Create 0 sized mesh to do all the generation of zero sized
624  // fields on processors that have zero sized meshes. Note that this is
625  // only necessary on master but since polyMesh construction with
626  // Pstream::parRun does parallel comms we have to do it on all
627  // processors
628  autoPtr<fvMeshSubset> subsetterPtr;
629 
630  // Missing a volume mesh somewhere?
631  if (volMeshOnProc.found(false))
632  {
633  // A zero-sized mesh with boundaries.
634  // This is used to create zero-sized fields.
635  subsetterPtr.reset(new fvMeshSubset(mesh, zero{}));
636  }
637 
638 
639  // Get original objects (before incrementing time!)
640  if (Pstream::master() && decompose)
641  {
642  runTime.caseName() = baseRunTime.caseName();
643  runTime.processorCase(false);
644  }
645  IOobjectList objects(mesh, runTime.timeName());
646  if (Pstream::master() && decompose)
647  {
648  runTime.caseName() = proc0CaseName;
649  runTime.processorCase(oldProcCase);
650  }
651 
652  Info<< "From time " << runTime.timeName()
653  << " mesh:" << mesh.objectRegistry::objectRelPath()
654  << " have objects:" << objects.names() << endl;
655 
656  // We don't want to map the decomposition (mapping already tested when
657  // mapping the cell centre field)
658  auto iter = objects.find("cellDist");
659  if (iter.good())
660  {
661  objects.erase(iter);
662  }
663 
664 
665  if (Pstream::master() && decompose)
666  {
667  runTime.caseName() = baseRunTime.caseName();
668  runTime.processorCase(false);
669  }
670 
671  // Field reading
672 
673  #undef doFieldReading
674  #define doFieldReading(Storage) \
675  { \
676  fieldsDistributor::readFields \
677  ( \
678  volMeshOnProc, mesh, subsetterPtr, objects, Storage \
679  ); \
680  }
681 
682  // volField
683  doFieldReading(volScalarFields);
684  doFieldReading(volVectorFields);
685  doFieldReading(volSphereTensorFields);
686  doFieldReading(volSymmTensorFields);
687  doFieldReading(volTensorFields);
688 
689  // surfaceField
690  doFieldReading(surfScalarFields);
691  doFieldReading(surfVectorFields);
692  doFieldReading(surfSphereTensorFields);
693  doFieldReading(surfSymmTensorFields);
694  doFieldReading(surfTensorFields);
695 
696  // Dimensioned internal fields
697  doFieldReading(dimScalarFields);
698  doFieldReading(dimVectorFields);
699  doFieldReading(dimSphereTensorFields);
700  doFieldReading(dimSymmTensorFields);
701  doFieldReading(dimTensorFields);
702 
703  // pointFields
704  nPointFields = 0;
705 
706  #undef doFieldReading
707  #define doFieldReading(Storage) \
708  { \
709  fieldsDistributor::readFields \
710  ( \
711  volMeshOnProc, oldPointMesh, subsetterPtr, objects, Storage, \
712  true /* (deregister field) */ \
713  ); \
714  nPointFields += Storage.size(); \
715  }
716 
717  doFieldReading(pointScalarFields);
718  doFieldReading(pointVectorFields);
719  doFieldReading(pointSphTensorFields);
720  doFieldReading(pointSymmTensorFields);
721  doFieldReading(pointTensorFields);
722  #undef doFieldReading
723 
724 
725  // Done reading
726 
727  if (Pstream::master() && decompose)
728  {
729  runTime.caseName() = proc0CaseName;
730  runTime.processorCase(oldProcCase);
731  }
732  }
733 
734  // Save pointMesh information before any topology changes occur!
735  if (nPointFields)
736  {
737  pointDistributor.saveMeshPoints();
738  }
739 
740 
741  // Mesh distribution engine
742  fvMeshDistribute distributor(mesh);
743 
744  // Do all the distribution of mesh and fields
745  autoPtr<mapDistributePolyMesh> distMap = distributor.distribute(decomp);
746 
747  // Print some statistics
748  Info<< "After distribution:" << endl;
749  printMeshData(mesh);
750 
751  // Get other side of processor boundaries
752  do
753  {
754  #undef doCorrectCoupled
755  #define doCorrectCoupled(FieldType) \
756  correctCoupledBoundaryConditions<processorFvPatch, FieldType>(mesh);
757 
758  doCorrectCoupled(volScalarField);
759  doCorrectCoupled(volVectorField);
760  doCorrectCoupled(volSphericalTensorField);
761  doCorrectCoupled(volSymmTensorField);
762  doCorrectCoupled(volTensorField);
763  #undef doCorrectCoupled
764  }
765  while (false);
766 
767  // No update surface fields
768 
769 
770  // Map pointFields
771  if (nPointFields)
772  {
773  // Construct new pointMesh from distributed mesh
774  const pointMesh& newPointMesh = pointMesh::New(mesh);
775 
776  pointDistributor.resetTarget(newPointMesh, distMap());
777 
778  pointDistributor.distributeAndStore(pointScalarFields);
779  pointDistributor.distributeAndStore(pointVectorFields);
780  pointDistributor.distributeAndStore(pointSphTensorFields);
781  pointDistributor.distributeAndStore(pointSymmTensorFields);
782  pointDistributor.distributeAndStore(pointTensorFields);
783  }
784 
785 
786  // Set the minimum write precision
788 
789 
790  if (!overwrite)
791  {
792  ++runTime;
794  }
795  else
796  {
797  mesh.setInstance(volMeshInstance);
798  }
799 
800 
801  // Register mapDistributePolyMesh for automatic writing...
802  IOmapDistributePolyMeshRef distMapRef
803  (
804  IOobject
805  (
806  "procAddressing",
809  mesh.thisDb(),
812  ),
813  distMap()
814  );
815 
816 
817  if (reconstruct)
818  {
819  if (Pstream::master())
820  {
821  Info<< "Setting caseName to " << baseRunTime.caseName()
822  << " to write reconstructed mesh (and fields)." << endl;
823  runTime.caseName() = baseRunTime.caseName();
824  const bool oldProcCase(runTime.processorCase(false));
825 
826  mesh.write();
828 
829  // Now we've written all. Reset caseName on master
830  Info<< "Restoring caseName" << endl;
831  runTime.caseName() = proc0CaseName;
832  runTime.processorCase(oldProcCase);
833  }
834  }
835  else
836  {
837  auto oldHandler = fileOperation::fileHandler(writeHandler);
838 
839  mesh.write();
840 
841  writeHandler = fileOperation::fileHandler(oldHandler);
842 
844  }
845  Info<< "Written redistributed mesh to "
846  << mesh.facesInstance() << nl << endl;
847 
848 
849  if (decompose || reconstruct)
850  {
851  // Decompose (1 -> N) or reconstruct (N -> 1)
852  // so {boundary,cell,face,point}ProcAddressing have meaning
854  (
855  mesh,
856  distMap(),
857  decompose,
858  writeHandler
859  );
860  }
861  else
862  {
863  // Redistribute (N -> M)
864  // {boundary,cell,face,point}ProcAddressing would be incorrect
865  // - can either remove or redistribute previous
867  }
868 
869 
870  // Refinement data
871  {
872 
873  // Read refinement data
874  if (Pstream::master() && decompose)
875  {
876  runTime.caseName() = baseRunTime.caseName();
877  runTime.processorCase(false);
878  }
879  IOobject io
880  (
881  "dummy",
884  mesh,
888  );
889 
890  hexRef8Data refData(io);
891  if (Pstream::master() && decompose)
892  {
893  runTime.caseName() = proc0CaseName;
894  runTime.processorCase(oldProcCase);
895  }
896 
897  // Make sure all processors have valid data (since only some will
898  // read)
899  refData.sync(io);
900 
901  // Distribute
902  refData.distribute(distMap());
903 
904 
905  // Now we've read refinement data we can remove it
907 
908  if (reconstruct)
909  {
910  if (Pstream::master())
911  {
912  const bool oldParRun = Pstream::parRun(false);
913 
914  Info<< "Setting caseName to " << baseRunTime.caseName()
915  << " to write reconstructed refinement data." << endl;
916  runTime.caseName() = baseRunTime.caseName();
917  const bool oldProcCase(runTime.processorCase(false));
918 
919  refData.write();
920 
921  // Now we've written all. Reset caseName on master
922  Info<< "Restoring caseName" << endl;
923  runTime.caseName() = proc0CaseName;
924  runTime.processorCase(oldProcCase);
925 
926  Pstream::parRun(oldParRun);
927  }
928  }
929  else
930  {
931  refData.write();
932  }
933  }
934 
936  //{
937  // // Read sets
938  // if (Pstream::master() && decompose)
939  // {
940  // runTime.caseName() = baseRunTime.caseName();
941  // runTime.processorCase(false);
942  // }
943  // IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
944  //
945  // PtrList<cellSet> cellSets;
946  // ReadFields(objects, cellSets);
947  //
948  // if (Pstream::master() && decompose)
949  // {
950  // runTime.caseName() = proc0CaseName;
951  // runTime.processorCase(oldProcCase);
952  // }
953  //
954  // forAll(cellSets, i)
955  // {
956  // cellSets[i].distribute(distMap());
957  // }
958  //
959  // if (reconstruct)
960  // {
961  // if (Pstream::master())
962  // {
963  // Info<< "Setting caseName to " << baseRunTime.caseName()
964  // << " to write reconstructed refinement data." << endl;
965  // runTime.caseName() = baseRunTime.caseName();
966  // const bool oldProcCase(runTime.processorCase(false));
967  //
968  // forAll(cellSets, i)
969  // {
970  // cellSets[i].distribute(distMap());
971  // }
972  //
973  // // Now we've written all. Reset caseName on master
974  // Info<< "Restoring caseName" << endl;
975  // runTime.caseName() = proc0CaseName;
976  // runTime.processorCase(oldProcCase);
977  // }
978  // }
979  // else
980  // {
981  // forAll(cellSets, i)
982  // {
983  // cellSets[i].distribute(distMap());
984  // }
985  // }
986  //}
987 
988 
989  return distMap;
990 }
991 
992 
993 /*---------------------------------------------------------------------------*\
994  main
995 \*---------------------------------------------------------------------------*/
996 
997 int main(int argc, char *argv[])
998 {
1000  (
1001  "Redistribute decomposed mesh and fields according"
1002  " to the decomposeParDict settings.\n"
1003  "Optionally run in decompose/reconstruct mode"
1004  );
1005 
1006  argList::noFunctionObjects(); // Never use function objects
1007 
1008  // enable -constant ... if someone really wants it
1009  // enable -zeroTime to prevent accidentally trashing the initial fields
1010  timeSelector::addOptions(true, true);
1011 
1012  #include "addAllRegionOptions.H"
1013 
1014  #include "addOverwriteOption.H"
1015  argList::addBoolOption("decompose", "Decompose case");
1016  argList::addBoolOption("reconstruct", "Reconstruct case");
1019  (
1020  "Test without writing the decomposition. "
1021  "Changes -cellDist to only write volScalarField."
1022  );
1023  argList::addVerboseOption("Additional verbosity");
1025  (
1026  "cellDist",
1027  "Write cell distribution as a labelList - for use with 'manual' "
1028  "decomposition method or as a volScalarField for post-processing."
1029  );
1031  (
1032  "newTimes",
1033  "Only reconstruct new times (i.e. that do not exist already)"
1034  );
1036  (
1037  "Additional verbosity. (Can be used multiple times)"
1038  );
1040  (
1041  "no-finite-area",
1042  "Suppress finiteArea mesh/field handling",
1043  true // Advanced option
1044  );
1045 
1046  // Handle arguments
1047  // ~~~~~~~~~~~~~~~~
1048  // (replacement for setRootCase that does not abort)
1049 
1050  argList args(argc, argv);
1051 
1052 
1053  // As much as possible avoid synchronised operation. To be looked at more
1054  // closely for the three scenarios:
1055  // - decompose - reads on master (and from parent directory) and sends
1056  // dictionary to slaves
1057  // - distribute - reads on potentially a different number of processors
1058  // than it writes to
1059  // - reconstruct - reads parallel, write on master only and to parent
1060  // directory
1061  refPtr<fileOperation> writeHandler;
1062  if
1063  (
1064  fileHandler().type()
1065  != fileOperations::uncollatedFileOperation::typeName
1066  )
1067  {
1068  // Install 'uncollated' as fileHandler. Save old one in writeHandler.
1069  writeHandler =
1071  }
1072 
1073  // Switch off parallel synchronisation of cached time directories
1074  fileHandler().distributed(true);
1075 
1076  // File handler to be used for writing
1077  const fileOperation& fh
1078  (
1079  writeHandler
1080  ? writeHandler()
1081  : fileHandler()
1082  );
1083 
1084  // Make sure to call findTimes on all processors to force caching of
1085  // time directories
1086  (void)fileHandler().findTimes(args.path(), "constant");
1087 
1088  // Need this line since we don't include "setRootCase.H"
1089  #include "foamDlOpenLibs.H"
1090 
1091  const bool reconstruct = args.found("reconstruct");
1092  const bool writeCellDist = args.found("cellDist");
1093  const bool dryrun = args.dryRun();
1094  const bool newTimes = args.found("newTimes");
1095  const int optVerbose = args.verbose();
1096 
1097  const bool doFiniteArea = !args.found("no-finite-area");
1098  bool decompose = args.found("decompose");
1099  bool overwrite = args.found("overwrite");
1100 
1101  if (optVerbose)
1102  {
1103  // Report on output
1106  }
1107 
1108  // Disable NaN setting and floating point error trapping. This is to avoid
1109  // any issues inside the field redistribution inside fvMeshDistribute
1110  // which temporarily moves processor faces into existing patches. These
1111  // will now not have correct values. After all bits have been assembled
1112  // the processor fields will be restored but until then there might
1113  // be uninitialised values which might upset some patch field constructors.
1114  // Workaround by disabling floating point error trapping. TBD: have
1115  // actual field redistribution instead of subsetting inside
1116  // fvMeshDistribute.
1117  Foam::sigFpe::unset(true);
1118 
1119  const wordRes selectedFields;
1120  const wordRes selectedLagrangianFields;
1121 
1122 
1123  if (decompose)
1124  {
1125  Info<< "Decomposing case (like decomposePar)" << nl << endl;
1126  if (reconstruct)
1127  {
1129  << "Cannot specify both -decompose and -reconstruct"
1130  << exit(FatalError);
1131  }
1132  }
1133  else if (reconstruct)
1134  {
1135  Info<< "Reconstructing case (like reconstructParMesh)" << nl << endl;
1136  }
1137 
1138  if ((decompose || reconstruct) && !overwrite)
1139  {
1140  overwrite = true;
1142  << "Working in -decompose or -reconstruct mode:"
1143  " automatically implies -overwrite" << nl << endl;
1144  }
1145 
1146  if (!Pstream::parRun())
1147  {
1149  << ": This utility can only be run parallel"
1150  << exit(FatalError);
1151  }
1152 
1153 
1154  if (!Foam::isDir(args.rootPath()))
1155  {
1157  << ": cannot open root directory " << args.rootPath()
1158  << exit(FatalError);
1159  }
1160 
1161  // Detect if running data-distributed (multiple roots)
1162  bool nfs = true;
1163  {
1164  List<fileName> roots(1, args.rootPath());
1166  nfs = (roots.size() == 1);
1167  }
1168 
1169  if (!nfs)
1170  {
1171  Info<< "Detected multiple roots i.e. non-nfs running"
1172  << nl << endl;
1173  }
1174 
1175  // Check if we have processor directories. Ideally would like to
1176  // use fileHandler().dirPath here but we don't have runTime yet and
1177  // want to delay constructing runTime until we've synced all time
1178  // directories...
1179  const fileName procDir(fileHandler().filePath(args.path()));
1180  if (Foam::isDir(procDir))
1181  {
1182  if (decompose)
1183  {
1184  Info<< "Removing existing processor directory:"
1185  << args.relativePath(procDir) << endl;
1186  fileHandler().rmDir(procDir);
1187  }
1188  }
1189  else
1190  {
1191  // Directory does not exist. If this happens on master -> decompose mode
1192  if (Pstream::master())
1193  {
1194  decompose = true;
1195  Info<< "No processor directories; switching on decompose mode"
1196  << nl << endl;
1197  }
1198  }
1199  // If master changed to decompose mode make sure all nodes know about it
1200  Pstream::broadcast(decompose);
1201 
1202 
1203  // If running distributed we have problem of new processors not finding
1204  // a system/controlDict. However if we switch on the master-only reading
1205  // the problem becomes that the time directories are differing sizes and
1206  // e.g. latestTime will pick up a different time (which causes createTime.H
1207  // to abort). So for now make sure to have master times on all
1208  // processors
1209  if (!reconstruct)
1210  {
1211  Info<< "Creating time directories on all processors" << nl << endl;
1212  createTimeDirs(args.path());
1213  }
1214 
1215  // Construct time
1216  // ~~~~~~~~~~~~~~
1217 
1218  #include "createTime.H"
1219  runTime.functionObjects().off(); // Extra safety?
1220 
1221 
1222  // Save local processor0 casename
1223  const fileName proc0CaseName = runTime.caseName();
1224  const bool oldProcCase = runTime.processorCase();
1225 
1226 
1227  // Construct undecomposed Time
1228  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1229  // This will read the same controlDict but might have a different
1230  // set of times so enforce same times
1231 
1232  if (!nfs)
1233  {
1234  Info<< "Creating time directories for undecomposed Time"
1235  << " on all processors" << nl << endl;
1236  createTimeDirs(args.globalPath());
1237  }
1238 
1239 
1240  Info<< "Create undecomposed database" << nl << endl;
1241  Time baseRunTime
1242  (
1243  runTime.controlDict(),
1244  runTime.rootPath(),
1246  runTime.system(),
1247  runTime.constant(),
1248  false, // No function objects
1249  false // No extra controlDict libs
1250  );
1251 
1252 
1253  wordHashSet masterTimeDirSet;
1254  if (newTimes)
1255  {
1256  instantList baseTimeDirs(baseRunTime.times());
1257  for (const instant& t : baseTimeDirs)
1258  {
1259  masterTimeDirSet.insert(t.name());
1260  }
1261  }
1262 
1263 
1264  // Allow override of decomposeParDict location
1265  const fileName decompDictFile =
1266  args.getOrDefault<fileName>("decomposeParDict", "");
1267 
1268  // Get region names
1269  #include "getAllRegionOptions.H"
1270 
1272  {
1273  Info<< "Using region: " << regionNames[0] << nl << endl;
1274  }
1275 
1276 
1277  // Demand driven lagrangian mapper
1278  autoPtr<parLagrangianDistributor> lagrangianDistributorPtr;
1279 
1280 
1281  if (reconstruct)
1282  {
1283  // use the times list from the master processor
1284  // and select a subset based on the command-line options
1286  Pstream::broadcast(timeDirs);
1287 
1288  if (timeDirs.empty())
1289  {
1291  << "No times selected"
1292  << exit(FatalError);
1293  }
1294 
1295 
1296  // Pass1 : reconstruct mesh and addressing
1297  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1298 
1299 
1300  Info<< nl
1301  << "Reconstructing mesh and addressing" << nl << endl;
1302 
1303  forAll(regionNames, regioni)
1304  {
1305  const word& regionName = regionNames[regioni];
1307 
1308  const fileName volMeshSubDir(regionDir/polyMesh::meshSubDir);
1309  const fileName areaMeshSubDir(regionDir/faMesh::meshSubDir);
1310 
1311  Info<< nl
1312  << "Reconstructing mesh " << regionDir << nl << endl;
1313 
1314  bool areaMeshDetected = false;
1315 
1316  // Loop over all times
1317  forAll(timeDirs, timeI)
1318  {
1319  // Set time for global database
1320  runTime.setTime(timeDirs[timeI], timeI);
1321  baseRunTime.setTime(timeDirs[timeI], timeI);
1322 
1323  Info<< "Time = " << runTime.timeName() << endl << endl;
1324 
1325  // Where meshes are
1326  fileName volMeshInstance;
1327  fileName areaMeshInstance;
1328 
1329  volMeshInstance = runTime.findInstance
1330  (
1331  volMeshSubDir,
1332  "faces",
1334  );
1335 
1336  if (doFiniteArea)
1337  {
1338  areaMeshInstance = runTime.findInstance
1339  (
1340  areaMeshSubDir,
1341  "faceLabels",
1343  );
1344  }
1345 
1347  (
1349  volMeshInstance,
1350  areaMeshInstance
1351  );
1352 
1353 
1354  // Check processors have meshes
1355  // - check for 'faces' file (polyMesh)
1356  // - check for 'faceLabels' file (faMesh)
1357  boolList volMeshOnProc;
1358  boolList areaMeshOnProc;
1359 
1360  volMeshOnProc = haveMeshFile
1361  (
1362  runTime,
1363  volMeshInstance/volMeshSubDir,
1364  "faces"
1365  );
1366 
1367  if (doFiniteArea)
1368  {
1369  areaMeshOnProc = haveMeshFile
1370  (
1371  runTime,
1372  areaMeshInstance/areaMeshSubDir,
1373  "faceLabels"
1374  );
1375 
1376  areaMeshDetected = areaMeshOnProc.found(true);
1377  }
1378 
1379 
1380  // Addressing back to reconstructed mesh as xxxProcAddressing.
1381  // - all processors have consistent faceProcAddressing
1382  // - processors without a mesh don't need faceProcAddressing
1383 
1384 
1385  // Note: filePath searches up on processors that don't have
1386  // processor if instance = constant so explicitly check
1387  // found filename.
1388  bool haveVolAddressing = false;
1389  if (volMeshOnProc[Pstream::myProcNo()])
1390  {
1391  // Read faces (just to know their size)
1392  faceCompactIOList faces
1393  (
1394  IOobject
1395  (
1396  "faces",
1397  volMeshInstance,
1398  volMeshSubDir,
1399  runTime,
1401  )
1402  );
1403 
1404  // Check faceProcAddressing
1405  labelIOList faceProcAddressing
1406  (
1407  IOobject
1408  (
1409  "faceProcAddressing",
1410  volMeshInstance,
1411  volMeshSubDir,
1412  runTime,
1414  )
1415  );
1416 
1417  haveVolAddressing =
1418  (
1419  faceProcAddressing.headerOk()
1420  && faceProcAddressing.size() == faces.size()
1421  );
1422  }
1423  else
1424  {
1425  // Have no mesh. Don't need addressing
1426  haveVolAddressing = true;
1427  }
1428 
1429  bool haveAreaAddressing = false;
1430  if (areaMeshOnProc[Pstream::myProcNo()])
1431  {
1432  // Read faces (just to know their size)
1434  (
1435  IOobject
1436  (
1437  "faceLabels",
1438  areaMeshInstance,
1439  areaMeshSubDir,
1440  runTime,
1442  )
1443  );
1444 
1445  // Check faceProcAddressing
1446  labelIOList faceProcAddressing
1447  (
1448  IOobject
1449  (
1450  "faceProcAddressing",
1451  areaMeshInstance,
1452  areaMeshSubDir,
1453  runTime,
1455  )
1456  );
1457 
1458  haveAreaAddressing =
1459  (
1460  faceProcAddressing.headerOk()
1461  && faceProcAddressing.size() == faceLabels.size()
1462  );
1463  }
1464  else if (areaMeshDetected)
1465  {
1466  // Have no mesh. Don't need addressing
1467  haveAreaAddressing = true;
1468  }
1469 
1470 
1471  // Additionally check for master faces being readable. Could
1472  // do even more checks, e.g. global number of cells same
1473  // as cellProcAddressing
1474 
1475  bool volMeshHaveUndecomposed = false;
1476  bool areaMeshHaveUndecomposed = false;
1477 
1478  if (Pstream::master())
1479  {
1480  Info<< "Checking " << baseRunTime.caseName()
1481  << " for undecomposed volume and area meshes..."
1482  << endl;
1483 
1484  const bool oldParRun = Pstream::parRun(false);
1485 
1486  // Volume
1487  {
1488  faceCompactIOList facesIO
1489  (
1490  IOobject
1491  (
1492  "faces",
1493  volMeshInstance,
1494  volMeshSubDir,
1495  baseRunTime,
1497  ),
1498  label(0)
1499  );
1500  volMeshHaveUndecomposed = facesIO.headerOk();
1501  }
1502 
1503  // Area
1504  if (doFiniteArea)
1505  {
1506  labelIOList labelsIO
1507  (
1508  IOobject
1509  (
1510  "faceLabels",
1511  areaMeshInstance,
1512  areaMeshSubDir,
1513  baseRunTime,
1515  )
1516  );
1517  areaMeshHaveUndecomposed = labelsIO.headerOk();
1518  }
1519 
1520  Pstream::parRun(oldParRun); // Restore parallel state
1521  }
1522 
1524  (
1526  volMeshHaveUndecomposed,
1527  areaMeshHaveUndecomposed
1528  );
1529 
1530  // Report
1531  {
1532  Info<< " volume mesh ["
1533  << volMeshHaveUndecomposed << "] : "
1534  << volMeshInstance << nl
1535  << " area mesh ["
1536  << areaMeshHaveUndecomposed << "] : "
1537  << areaMeshInstance << nl
1538  << endl;
1539  }
1540 
1541 
1542  if
1543  (
1544  !volMeshHaveUndecomposed
1545  || !returnReduceAnd(haveVolAddressing)
1546  )
1547  {
1548  Info<< "No undecomposed mesh. Creating from: "
1549  << volMeshInstance << endl;
1550 
1551  if (areaMeshHaveUndecomposed)
1552  {
1553  areaMeshHaveUndecomposed = false;
1554  Info<< "Also ignore any undecomposed area mesh"
1555  << endl;
1556  }
1557 
1559  (
1560  IOobject
1561  (
1562  regionName,
1563  volMeshInstance,
1564  runTime,
1566  ),
1567  decompose
1568  );
1569  fvMeshTools::setBasicGeometry(volMeshPtr());
1570  fvMesh& mesh = volMeshPtr();
1571 
1572 
1573  Info<< nl << "Reconstructing mesh" << nl << endl;
1574 
1575  // Reconstruct (1 processor)
1576  const label nDestProcs(1);
1577  const labelList finalDecomp(mesh.nCells(), Zero);
1578 
1579  redistributeAndWrite
1580  (
1581  writeHandler,
1582  baseRunTime,
1583  proc0CaseName,
1584 
1585  // Controls
1586  false, // do not read fields
1587  false, // do not read undecomposed case on proc0
1588  true, // write redistributed files to proc0
1589  overwrite,
1590 
1591  // Decomposition information
1592  nDestProcs,
1593  finalDecomp,
1594 
1595  // For finite-volume
1596  volMeshOnProc,
1597  volMeshInstance,
1598  mesh
1599  );
1600  }
1601 
1602 
1603  // Similarly for finiteArea
1604  // - may or may not have undecomposed mesh
1605  // - may or may not have decomposed meshes
1606 
1607  if
1608  (
1609  areaMeshOnProc.found(true) // ie, areaMeshDetected
1610  &&
1611  (
1612  !areaMeshHaveUndecomposed
1613  || !returnReduceAnd(haveAreaAddressing)
1614  )
1615  )
1616  {
1617  Info<< "Loading area mesh from "
1618  << areaMeshInstance << endl;
1619 
1620  Info<< " getting volume mesh support" << endl;
1621 
1623  (
1624  IOobject
1625  (
1626  regionName,
1627  baseRunTime.timeName(),
1628  baseRunTime,
1630  ),
1631  true // read on master only
1632  );
1633  fvMeshTools::setBasicGeometry(baseMeshPtr());
1634 
1636  (
1637  IOobject
1638  (
1639  regionName,
1640  baseMeshPtr().facesInstance(),
1641  runTime,
1643  ),
1644  decompose
1645  );
1646  fvMeshTools::setBasicGeometry(volMeshPtr());
1647  fvMesh& mesh = volMeshPtr();
1648 
1649  // Read volume proc addressing back to base mesh
1651  (
1653  );
1654 
1655 
1657  (
1658  IOobject
1659  (
1660  regionName,
1661  areaMeshInstance,
1662  runTime,
1664  ),
1665  mesh, // <- The referenced polyMesh (from above)
1666  decompose
1667  );
1668  faMesh& areaMesh = areaMeshPtr();
1669 
1672 
1673  autoPtr<faMesh> areaBaseMeshPtr;
1674 
1675  // Reconstruct using polyMesh distribute map
1676  mapDistributePolyMesh faDistMap
1677  (
1679  (
1680  areaMesh,
1681  distMap(), // The polyMesh distMap
1682  baseMeshPtr(), // Target polyMesh
1683  areaBaseMeshPtr
1684  )
1685  );
1686 
1687  faMeshTools::forceDemandDriven(areaBaseMeshPtr());
1688  faMeshTools::unregisterMesh(areaBaseMeshPtr());
1689 
1690 
1691  if (Pstream::master())
1692  {
1693  Info<< "Setting caseName to " << baseRunTime.caseName()
1694  << " to write reconstructed area mesh." << endl;
1695  runTime.caseName() = baseRunTime.caseName();
1696  const bool oldProcCase(runTime.processorCase(false));
1697 
1698  areaBaseMeshPtr().write();
1699 
1700  // Now we've written all. Reset caseName on master
1701  Info<< "Restoring caseName" << endl;
1702  runTime.caseName() = proc0CaseName;
1703  runTime.processorCase(oldProcCase);
1704  }
1705 
1706  // Update for the reconstructed procAddressing
1708  (
1709  areaBaseMeshPtr(), // Reconstruct location
1710  faDistMap,
1711  false, // decompose=false
1712  writeHandler,
1713  areaMeshPtr.get() // procMesh
1714  );
1715  }
1716  }
1717 
1718  // Make sure all is finished writing until re-reading in pass2
1719  // below
1720  fileHandler().flush();
1721 
1722 
1723  // Pass2 : read mesh and addressing and reconstruct fields
1724  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1725 
1726  Info<< nl
1727  << "Reconstructing fields" << nl << endl;
1728 
1729  runTime.setTime(timeDirs[0], 0);
1730  baseRunTime.setTime(timeDirs[0], 0);
1731  Info<< "Time = " << runTime.timeName() << endl << endl;
1732 
1733 
1734  // Read undecomposed mesh on master and 'empty' mesh
1735  // (zero faces, point, cells but valid patches and zones) on slaves.
1736  // This is a bit of tricky code and hidden inside fvMeshTools for
1737  // now.
1738  Info<< "Reading undecomposed mesh (on master)" << endl;
1739 
1741  (
1742  IOobject
1743  (
1744  regionName,
1745  baseRunTime.timeName(),
1746  baseRunTime,
1748  ),
1749  true // read on master only
1750  );
1751 
1752  fvMeshTools::setBasicGeometry(baseMeshPtr());
1753 
1754  Info<< "Reading local, decomposed mesh" << endl;
1756  (
1757  IOobject
1758  (
1759  regionName,
1760  baseMeshPtr().facesInstance(),
1761  runTime,
1763  ),
1764  decompose
1765  );
1766  fvMesh& mesh = volMeshPtr();
1767 
1768 
1769  // Similarly for finiteArea
1770  autoPtr<faMesh> areaBaseMeshPtr;
1771  autoPtr<faMesh> areaMeshPtr;
1772  autoPtr<faMeshDistributor> faDistributor;
1773  mapDistributePolyMesh areaDistMap;
1774 
1775  if (areaMeshDetected)
1776  {
1777  areaBaseMeshPtr = faMeshTools::newMesh
1778  (
1779  IOobject
1780  (
1781  regionName,
1782  baseRunTime.timeName(),
1783  baseRunTime,
1785  ),
1786  baseMeshPtr(),
1787  true // read on master only
1788  );
1789 
1790  areaMeshPtr = faMeshTools::loadOrCreateMesh
1791  (
1792  IOobject
1793  (
1794  regionName,
1795  areaBaseMeshPtr().facesInstance(),
1796  runTime,
1798  ),
1799  mesh,
1800  decompose
1801  );
1802 
1803  areaDistMap =
1805  (
1806  areaMeshPtr(),
1807  areaBaseMeshPtr
1808  );
1809 
1810  faMeshTools::forceDemandDriven(areaMeshPtr());
1811 
1812  // Create an appropriate field distributor
1813  faDistributor.reset
1814  (
1815  new faMeshDistributor
1816  (
1817  areaMeshPtr(), // source
1818  areaBaseMeshPtr(), // target
1819  areaDistMap,
1820  Pstream::master() // only write on master
1821  )
1822  );
1823  // Report some messages. Tbd.
1825  }
1826 
1827 
1828  if (writeHandler && Pstream::master())
1829  {
1830  // Remove any left-over empty processor directories created
1831  // by loadOrCreateMesh to get around e.g. collated start-up
1832  // problems. Should not happen in reconstruct mode ...
1833  const bool oldParRun = Pstream::parRun(false);
1835  Pstream::parRun(oldParRun);
1836  }
1837 
1838 
1839  // Read addressing back to base mesh
1841  distMap = fvMeshTools::readProcAddressing(mesh, baseMeshPtr);
1842 
1843  // Construct field mapper
1844  auto fvDistributorPtr =
1846  (
1847  mesh, // source
1848  baseMeshPtr(), // target
1849  distMap(),
1850  UPstream::master() // Write reconstructed on master
1851  );
1852 
1853  // Construct point field mapper
1854  const auto& basePointMesh = pointMesh::New(baseMeshPtr());
1855  const auto& procPointMesh = pointMesh::New(mesh);
1856 
1857  auto pointFieldDistributorPtr =
1859  (
1860  procPointMesh, // source
1861  basePointMesh, // target
1862  distMap(),
1863  false, // delay
1864  UPstream::master() // Write reconstructed on master
1865  );
1866 
1867 
1868  // Since we start from Times[0] and not runTime.timeName() we
1869  // might overlook point motion in the first timestep
1870  // (since mesh.readUpdate() below will not be triggered). Instead
1871  // detect points by hand
1873  {
1874  Info<< " Detected initial mesh motion;"
1875  << " reconstructing points" << nl
1876  << endl;
1877  fvDistributorPtr().reconstructPoints();
1878  }
1879 
1880 
1881  // Loop over all times
1882  forAll(timeDirs, timeI)
1883  {
1884  if (newTimes && masterTimeDirSet.found(timeDirs[timeI].name()))
1885  {
1886  Info<< "Skipping time " << timeDirs[timeI].name()
1887  << nl << endl;
1888  continue;
1889  }
1890 
1891  // Set time for global database
1892  runTime.setTime(timeDirs[timeI], timeI);
1893  baseRunTime.setTime(timeDirs[timeI], timeI);
1894 
1895  Info<< "Time = " << runTime.timeName() << endl << endl;
1896 
1897 
1898  // Check if any new meshes need to be read.
1900 
1901  if (procStat == polyMesh::POINTS_MOVED)
1902  {
1903  Info<< " Detected mesh motion; reconstructing points"
1904  << nl << endl;
1905  fvDistributorPtr().reconstructPoints();
1906  }
1907  else if
1908  (
1909  procStat == polyMesh::TOPO_CHANGE
1910  || procStat == polyMesh::TOPO_PATCH_CHANGE
1911  )
1912  {
1913  Info<< " Detected topology change;"
1914  << " reconstructing addressing" << nl << endl;
1915 
1916  if (baseMeshPtr)
1917  {
1918  // Cannot do a baseMesh::readUpdate() since not all
1919  // processors will have mesh files. So instead just
1920  // recreate baseMesh
1921  baseMeshPtr.clear();
1922  baseMeshPtr = fvMeshTools::newMesh
1923  (
1924  IOobject
1925  (
1926  regionName,
1927  baseRunTime.timeName(),
1928  baseRunTime,
1930  ),
1931  true // read on master only
1932  );
1933  }
1934 
1935  // Re-read procAddressing
1936  distMap =
1937  fvMeshTools::readProcAddressing(mesh, baseMeshPtr);
1938 
1939  // Reset field mappers
1940 
1941  fvDistributorPtr.reset
1942  (
1944  (
1945  mesh, // source
1946  baseMeshPtr(), // target
1947  distMap(),
1948  UPstream::master() // Write reconstruct on master
1949  )
1950  );
1951 
1952  // Construct point field mapper
1953  const auto& basePointMesh = pointMesh::New(baseMeshPtr());
1954  const auto& procPointMesh = pointMesh::New(mesh);
1955 
1956  pointFieldDistributorPtr.reset
1957  (
1959  (
1960  procPointMesh, // source
1961  basePointMesh, // target
1962  distMap(),
1963  false, // delay until later
1964  UPstream::master() // Write reconstruct on master
1965  )
1966  );
1967 
1968  lagrangianDistributorPtr.reset();
1969 
1970  if (areaMeshPtr)
1971  {
1972  Info<< " Discarding finite-area addressing"
1973  << " (TODO)" << nl << endl;
1974 
1975  areaBaseMeshPtr.reset();
1976  areaMeshPtr.reset();
1977  faDistributor.reset();
1978  areaDistMap.clear();
1979  }
1980  }
1981 
1982 
1983  // Get list of objects
1984  IOobjectList objects(mesh, runTime.timeName());
1985 
1986 
1987  // Mesh fields (vol, surface, volInternal)
1988  fvDistributorPtr()
1989  .distributeAllFields(objects, selectedFields);
1990 
1991  // pointfields
1992  // - distribute and write (verbose)
1993  pointFieldDistributorPtr()
1994  .distributeAllFields(objects, selectedFields);
1995 
1996 
1997  // Clouds (note: might not be present on all processors)
1999  (
2000  lagrangianDistributorPtr,
2001  baseMeshPtr(),
2002  mesh,
2003  distMap(),
2004  selectedLagrangianFields
2005  );
2006 
2007  if (faDistributor)
2008  {
2009  faDistributor()
2010  .distributeAllFields(objects, selectedFields);
2011  }
2012 
2013  // If there are any "uniform" directories copy them from
2014  // the master processor
2015  copyUniform
2016  (
2017  fh,
2018  decompose,
2019  reconstruct,
2020  mesh.time().timeName(),
2021  mesh,
2022  baseMeshPtr()
2023  );
2024  // Non-region specific. Note: should do outside region loop
2025  // but would then have to replicate the whole time loop ...
2026  copyUniform
2027  (
2028  fh,
2029  decompose,
2030  reconstruct,
2031  mesh.time().timeName(),
2032  mesh.time(), // runTime
2033  baseMeshPtr().time() // baseRunTime
2034  );
2035  }
2036  }
2037  }
2038  else
2039  {
2040  // decompose or redistribution mode.
2041  // decompose : master : read from parent dir
2042  // slave : dummy mesh
2043  // redistribute : all read mesh or dummy mesh
2044 
2045  // Time coming from processor0 (or undecomposed if no processor0)
2046  scalar masterTime;
2047  if (decompose)
2048  {
2049  // Use base time. This is to handle e.g. startTime = latestTime
2050  // which will not do anything if there are no processor directories
2051  masterTime = timeSelector::selectIfPresent
2052  (
2053  baseRunTime,
2054  args
2055  )[0].value();
2056  }
2057  else
2058  {
2059  masterTime = timeSelector::selectIfPresent
2060  (
2061  runTime,
2062  args
2063  )[0].value();
2064  }
2065  Pstream::broadcast(masterTime);
2066  Info<< "Setting time to that of master or undecomposed case : "
2067  << masterTime << endl;
2068  runTime.setTime(masterTime, 0);
2069  baseRunTime.setTime(masterTime, 0);
2070 
2071  // Save old time name (since might be incremented)
2072  const word oldTimeName(runTime.timeName());
2073 
2074  forAll(regionNames, regioni)
2075  {
2076  const word& regionName = regionNames[regioni];
2078 
2079  const fileName volMeshSubDir(regionDir/polyMesh::meshSubDir);
2080  const fileName areaMeshSubDir(regionDir/faMesh::meshSubDir);
2081 
2082  Info<< nl << nl
2083  << (decompose ? "Decomposing" : "Redistributing")
2084  << " mesh " << regionDir << nl << endl;
2085 
2086 
2087  // Get time instance directory
2088  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
2089  // At this point we should be able to read at least a mesh on
2090  // processor0. Note the changing of the processor0 casename to
2091  // enforce it to read/write from the undecomposed case
2092 
2093  fileName volMeshMasterInstance;
2094  fileName areaMeshMasterInstance;
2095 
2096  // Assume to be true
2097  bool volMeshHaveUndecomposed = true;
2098  bool areaMeshHaveUndecomposed = doFiniteArea;
2099 
2100  if (Pstream::master())
2101  {
2102  if (decompose)
2103  {
2104  Info<< "Checking undecomposed mesh in case: "
2105  << baseRunTime.caseName() << endl;
2106  runTime.caseName() = baseRunTime.caseName();
2107  runTime.processorCase(false);
2108  }
2109 
2110  const bool oldParRun = Pstream::parRun(false);
2111  volMeshMasterInstance = runTime.findInstance
2112  (
2113  volMeshSubDir,
2114  "faces",
2116  );
2117 
2118  if (doFiniteArea)
2119  {
2120  areaMeshMasterInstance = runTime.findInstance
2121  (
2122  areaMeshSubDir,
2123  "faceLabels",
2125  );
2126 
2127  // Note: findInstance returns "constant" even if not found,
2128  // so recheck now for a false positive.
2129 
2130  if ("constant" == areaMeshMasterInstance)
2131  {
2132  const boolList areaMeshOnProc
2133  (
2134  haveMeshFile
2135  (
2136  runTime,
2137  areaMeshMasterInstance/areaMeshSubDir,
2138  "faceLabels",
2139  false // verbose=false
2140  )
2141  );
2142 
2143  if (areaMeshOnProc.empty() || !areaMeshOnProc[0])
2144  {
2145  areaMeshHaveUndecomposed = false;
2146  }
2147  }
2148  }
2149 
2150  Pstream::parRun(oldParRun); // Restore parallel state
2151 
2152  if (decompose)
2153  {
2154  Info<< " volume mesh ["
2155  << volMeshHaveUndecomposed << "] : "
2156  << volMeshMasterInstance << nl
2157  << " area mesh ["
2158  << areaMeshHaveUndecomposed << "] : "
2159  << areaMeshMasterInstance << nl
2160  << nl << nl;
2161 
2162  // Restoring caseName
2163  runTime.caseName() = proc0CaseName;
2164  runTime.processorCase(oldProcCase);
2165  }
2166  }
2167 
2169  (
2171  volMeshHaveUndecomposed,
2172  areaMeshHaveUndecomposed,
2173  volMeshMasterInstance,
2174  areaMeshMasterInstance
2175  );
2176 
2177  // Check processors have meshes
2178  // - check for 'faces' file (polyMesh)
2179  // - check for 'faceLabels' file (faMesh)
2180  boolList volMeshOnProc;
2181  boolList areaMeshOnProc;
2182 
2183  volMeshOnProc = haveMeshFile
2184  (
2185  runTime,
2186  volMeshMasterInstance/volMeshSubDir,
2187  "faces"
2188  );
2189 
2190  if (doFiniteArea)
2191  {
2192  areaMeshOnProc = haveMeshFile
2193  (
2194  runTime,
2195  areaMeshMasterInstance/areaMeshSubDir,
2196  "faceLabels"
2197  );
2198  }
2199 
2200  // Prior to loadOrCreateMesh, note which meshes already exist
2201  // for the current file handler.
2202  // - where mesh would be written if it didn't exist already.
2203  fileNameList volMeshDir(Pstream::nProcs());
2204  {
2205  volMeshDir[Pstream::myProcNo()] =
2206  (
2207  fileHandler().objectPath
2208  (
2209  IOobject
2210  (
2211  "faces",
2212  volMeshMasterInstance/volMeshSubDir,
2213  runTime
2214  ),
2215  word::null
2216  ).path()
2217  );
2218 
2219  Pstream::allGatherList(volMeshDir);
2220 
2221  if (optVerbose && Pstream::master())
2222  {
2223  Info<< "Per processor faces dirs:" << nl
2224  << '(' << nl;
2225 
2226  for (const int proci : Pstream::allProcs())
2227  {
2228  Info<< " "
2229  << runTime.relativePath(volMeshDir[proci]);
2230 
2231  if (!volMeshOnProc[proci])
2232  {
2233  Info<< " [missing]";
2234  }
2235  Info<< nl;
2236  }
2237  Info<< ')' << nl << endl;
2238  }
2239  }
2240 
2241  fileNameList areaMeshDir(Pstream::nProcs());
2242  if (doFiniteArea)
2243  {
2244  areaMeshDir[Pstream::myProcNo()] =
2245  (
2246  fileHandler().objectPath
2247  (
2248  IOobject
2249  (
2250  "faceLabels",
2251  areaMeshMasterInstance/areaMeshSubDir,
2252  runTime
2253  ),
2254  word::null
2255  ).path()
2256  );
2257 
2258  Pstream::allGatherList(areaMeshDir);
2259 
2260  if (optVerbose && Pstream::master())
2261  {
2262  Info<< "Per processor faceLabels dirs:" << nl
2263  << '(' << nl;
2264 
2265  for (const int proci : Pstream::allProcs())
2266  {
2267  Info<< " "
2268  << runTime.relativePath(areaMeshDir[proci]);
2269 
2270  if (!areaMeshOnProc[proci])
2271  {
2272  Info<< " [missing]";
2273  }
2274  Info<< nl;
2275  }
2276  Info<< ')' << nl << endl;
2277  }
2278  }
2279 
2280 
2281  // Load mesh (or create dummy one)
2282  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2283 
2284  if (Pstream::master() && decompose)
2285  {
2286  Info<< "Setting caseName to " << baseRunTime.caseName()
2287  << " to read undecomposed mesh" << endl;
2288  runTime.caseName() = baseRunTime.caseName();
2289  runTime.processorCase(false);
2290  }
2291 
2292  // Volume mesh
2294  (
2295  IOobject
2296  (
2297  regionName,
2298  volMeshMasterInstance,
2299  runTime,
2301  ),
2302  decompose
2303  );
2304  fvMesh& mesh = volMeshPtr();
2305 
2306 
2307  // Area mesh
2308 
2309  autoPtr<faMesh> areaMeshPtr;
2310 
2311  // Decomposing: must have an undecomposed mesh
2312  // Redistributing: have any proc mesh
2313  if
2314  (
2315  doFiniteArea
2316  &&
2317  (
2318  decompose
2319  ? areaMeshHaveUndecomposed
2320  : areaMeshOnProc.found(true)
2321  )
2322  )
2323  {
2324  areaMeshPtr = faMeshTools::loadOrCreateMesh
2325  (
2326  IOobject
2327  (
2328  regionName,
2329  areaMeshMasterInstance,
2330  runTime,
2332  ),
2333  mesh, // <- The referenced polyMesh (from above)
2334  decompose
2335  );
2336 
2337  faMeshTools::forceDemandDriven(*areaMeshPtr);
2338  faMeshTools::unregisterMesh(*areaMeshPtr);
2339  }
2340 
2341 
2342  if (writeHandler)
2343  {
2344  // Remove any left-over empty processor directories created
2345  // by loadOrCreateMesh to get around the collated start-up
2346  // problems
2347  Info<< "Removing left-over empty processor directories" << nl;
2348 
2349  if (Pstream::master()) //fileHandler().comm()))
2350  {
2351  const auto myProci = UPstream::myProcNo(); //comm()
2352  const auto& procs = UPstream::procID
2353  (
2355  );
2356  const bool oldParRun = Pstream::parRun(false);
2357  for (const auto proci : procs)
2358  {
2359  if
2360  (
2361  !volMeshOnProc[proci]
2362  && volMeshDir[proci] != volMeshDir[myProci]
2363  )
2364  {
2365  Foam::rmDir(volMeshDir[proci], true); // silent
2366  }
2367 
2368  if
2369  (
2370  !areaMeshOnProc[proci]
2371  && areaMeshDir[proci] != areaMeshDir[myProci]
2372  )
2373  {
2374  Foam::rmDir(areaMeshDir[proci], true); // silent
2375  }
2376 
2377  // Remove empty processor directories
2378  // Eg, <path-name>/processorN/constant/polyMesh
2379  // to <path-name>/processorN
2380  if (proci != myProci)
2381  {
2383  (
2384  volMeshDir[proci].path().path()
2385  );
2386  }
2387  }
2388 
2389  // Remove empty directory
2391 
2392  Pstream::parRun(oldParRun);
2393  }
2394  }
2395 
2396 
2397  if (Pstream::master() && decompose)
2398  {
2399  Info<< "Restoring caseName" << endl;
2400  runTime.caseName() = proc0CaseName;
2401  runTime.processorCase(oldProcCase);
2402  }
2403 
2404  const label nOldCells = mesh.nCells();
2405 
2406  // const label nOldAreaFaces =
2407  // (areaMeshPtr ? areaMeshPtr().nFaces() : 0);
2408  //
2409  //Pout<< "Loaded mesh : nCells:" << nOldCells
2410  // << " nPatches:" << mesh.boundaryMesh().size() << endl;
2411  //Pout<< "Loaded area mesh : nFaces:" << nOldAreaFaces << endl;
2412 
2413  // Determine decomposition
2414  // ~~~~~~~~~~~~~~~~~~~~~~~
2415 
2416  label nDestProcs;
2417  labelList finalDecomp;
2418  determineDecomposition
2419  (
2420  baseRunTime,
2421  decompDictFile,
2422  decompose,
2423  proc0CaseName,
2424  mesh,
2425  writeCellDist,
2426 
2427  nDestProcs,
2428  finalDecomp
2429  );
2430 
2431 
2432  if (dryrun)
2433  {
2434  if (!Pstream::master())
2435  {
2436  if (areaMeshPtr && !areaMeshOnProc[Pstream::myProcNo()])
2437  {
2438  // Remove dummy mesh created by loadOrCreateMesh
2439  const bool oldParRun = Pstream::parRun(false);
2440  areaMeshPtr->removeFiles();
2441  Pstream::parRun(oldParRun); // Restore parallel state
2442  }
2443 
2444  if (!volMeshOnProc[Pstream::myProcNo()])
2445  {
2446  // Remove dummy mesh created by loadOrCreateMesh
2447  const bool oldParRun = Pstream::parRun(false);
2448  mesh.removeFiles();
2449  // Silent rmdir
2450  Foam::rmDir(mesh.objectRegistry::objectPath(), true);
2451  Pstream::parRun(oldParRun); // Restore parallel state
2452  }
2453  }
2454  continue;
2455  }
2456 
2457  // Area fields first. Read and deregister
2459  if (areaMeshPtr)
2460  {
2461  areaFields.read
2462  (
2463  baseRunTime,
2464  proc0CaseName,
2465  decompose,
2466 
2467  areaMeshOnProc,
2468  areaMeshMasterInstance,
2469  (*areaMeshPtr)
2470  );
2471  }
2472 
2473 
2474  // Detect lagrangian fields
2475  if (Pstream::master() && decompose)
2476  {
2477  runTime.caseName() = baseRunTime.caseName();
2478  runTime.processorCase(false);
2479  }
2480 
2481  // Read lagrangian fields and store on cloud (objectRegistry)
2483  (
2485  (
2486  mesh,
2487  selectedLagrangianFields
2488  )
2489  );
2490 
2491  if (Pstream::master() && decompose)
2492  {
2493  runTime.caseName() = proc0CaseName;
2494  runTime.processorCase(oldProcCase);
2495  }
2496 
2497  // Load fields, do all distribution (mesh and fields)
2498  // - but not lagrangian fields; these are done later
2499  autoPtr<mapDistributePolyMesh> distMap = redistributeAndWrite
2500  (
2501  writeHandler,
2502  baseRunTime,
2503  proc0CaseName,
2504 
2505  // Controls
2506  true, // read fields
2507  decompose, // decompose, i.e. read from undecomposed case
2508  false, // no reconstruction
2509  overwrite,
2510 
2511  // Decomposition information
2512  nDestProcs,
2513  finalDecomp,
2514 
2515  // For finite volume
2516  volMeshOnProc,
2517  volMeshMasterInstance,
2518  mesh
2519  );
2520 
2521 
2522  // Redistribute any clouds
2524  (
2525  lagrangianDistributorPtr,
2526  mesh,
2527  nOldCells,
2528  distMap(),
2529  clouds
2530  );
2531 
2532 
2533  // Redistribute area fields
2534 
2535  mapDistributePolyMesh faDistMap;
2536  autoPtr<faMesh> areaProcMeshPtr;
2537 
2538  if (areaMeshPtr)
2539  {
2540  faDistMap = faMeshDistributor::distribute
2541  (
2542  areaMeshPtr(),
2543  distMap(),
2544  areaProcMeshPtr
2545  );
2546 
2547  // Force recreation of everything that might vaguely
2548  // be used by patches:
2549 
2550  faMeshTools::forceDemandDriven(areaProcMeshPtr());
2551 
2552 
2553  if (reconstruct)
2554  {
2555  if (Pstream::master())
2556  {
2557  Info<< "Setting caseName to " << baseRunTime.caseName()
2558  << " to write reconstructed mesh (and fields)."
2559  << endl;
2560  runTime.caseName() = baseRunTime.caseName();
2561  const bool oldProcCase(runTime.processorCase(false));
2562  //const bool oldParRun = Pstream::parRun(false);
2563 
2564  areaProcMeshPtr->write();
2565 
2566  // Now we've written all. Reset caseName on master
2567  Info<< "Restoring caseName" << endl;
2568  runTime.caseName() = proc0CaseName;
2569  runTime.processorCase(oldProcCase);
2570  }
2571  }
2572  else
2573  {
2574  auto oldHandler = fileOperation::fileHandler(writeHandler);
2575 
2577  (
2578  IOobject
2579  (
2580  "procAddressing",
2581  areaProcMeshPtr->facesInstance(),
2583  areaProcMeshPtr->thisDb(),
2587  ),
2588  faDistMap
2589  ).write();
2590 
2591  areaProcMeshPtr->write();
2592 
2593  writeHandler = fileOperation::fileHandler(oldHandler);
2594 
2595  if (decompose)
2596  {
2598  (
2599  areaProcMeshPtr(),
2600  faDistMap,
2601  decompose,
2602  writeHandler
2603  );
2604  }
2605  }
2606 
2607  Info<< "Written redistributed mesh to "
2608  << areaProcMeshPtr->facesInstance() << nl << endl;
2609 
2610  faMeshDistributor distributor
2611  (
2612  areaMeshPtr(), // source
2613  areaProcMeshPtr(), // target
2614  faDistMap,
2615  true // isWriteProc (unused)
2616  );
2617 
2618  areaFields.redistributeAndWrite(distributor, true);
2619  }
2620 
2621  // Copy region-specific uniform
2622  // (e.g. solid/uniform/cumulativeContErr)
2623  copyUniform
2624  (
2625  fh,
2626  decompose,
2627  reconstruct,
2628  oldTimeName, // provided read time
2629  mesh, // read location is mesh (but oldTimeName)
2630  mesh // write location is mesh
2631  );
2632  }
2633 
2634  // Copy non-region specific uniform (e.g. uniform/time)
2635  copyUniform
2636  (
2637  fh,
2638  decompose,
2639  reconstruct,
2640  oldTimeName, // provided read time
2641  ( // read location
2642  decompose
2643  ? baseRunTime
2644  : runTime
2645  ),
2646  runTime // writing location
2647  );
2648  }
2649 
2650 
2651  Info<< "End\n" << endl;
2652 
2653  return 0;
2654 }
2655 
2656 
2657 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:218
Foam::surfaceFields.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:88
static autoPtr< faMesh > newMesh(const IOobject &io, const polyMesh &pMesh, const bool masterOnlyReading, const bool verbose=false)
Read mesh or create dummy mesh (0 faces, >0 patches).
Definition: faMeshTools.C:66
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
virtual fileName dirPath(const bool checkGlobal, const IOobject &io, const bool search=true) const =0
Search for a directory. checkGlobal : also check undecomposed.
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
boolList haveMeshFile(const Time &runTime, const fileName &meshPath, const word &meshFile="faces", const bool verbose=true)
Check for availability of specified mesh file (default: "faces")
fileName path() const
Return path.
Definition: Time.H:454
A class for handling file names.
Definition: fileName.H:71
Inter-processor communication reduction functions.
virtual bool parallelAware() const =0
Is method parallel aware?
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:854
void off()
Switch the function objects off.
word findInstance(const fileName &dir, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null) const
Return time instance (location) of dir that contains the file name (eg, used in reading mesh data)...
Definition: Time.C:787
static void writeProcAddressing(const faMesh &mesh, const mapDistributePolyMesh &faDistMap, const bool decompose, refPtr< fileOperation > &writeHandler, const faMesh *procMesh=nullptr)
Write decompose/reconstruct addressing.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
static void removeFiles(const polyMesh &)
Helper: remove all relevant files from mesh instance.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
static void unset(bool verbose=false)
Deactivate SIGFPE signal handler and NaN memory initialisation.
Definition: sigFpe.C:234
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:842
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:1131
const fileName & caseName() const noexcept
Return case name.
Definition: TimePathsI.H:55
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
int debugSwitch(const char *name, const int deflt=0)
Lookup debug switch or add default value.
Definition: debug.C:222
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:632
static label nDomains(const dictionary &decompDict, const word &regionName="")
Return region-specific or top-level numberOfSubdomains entry.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:402
wordList names() const
The unsorted names of all objects.
virtual labelList decompose(const pointField &points, const scalarField &pointWeights) const
Return the wanted processor number for every coordinate.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
bool found(const T &val, label pos=0) const
Same as contains()
Definition: UList.H:840
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1004
void clear() noexcept
Same as reset(nullptr)
Definition: autoPtr.H:255
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:80
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:423
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
const fileName & rootPath() const noexcept
Return root path.
Definition: TimePathsI.H:43
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
wordList regionNames
void clear()
Reset to zero size, only retaining communicator(s)
Ignore writing from objectRegistry::writeObject()
Reading, reconstruct, redistribution of lagrangian fields.
const dimensionSet dimless
Dimensionless.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1029
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
static const decompositionModel & New(const polyMesh &mesh, const fileName &decompDictFile="", const dictionary *fallback=nullptr)
Read and register on mesh, optionally with alternative decomposeParDict path/name or with fallback co...
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:362
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:411
A IOmapDistributePolyMesh wrapper for using referenced external data.
labelList faceLabels(nFaceLabels)
Distributor/redistributor for point fields, uses a two (or three) stage construction.
decompositionMethod & decomposer() const
Return demand-driven decomposition method.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
instantList select(const instantList &times) const
Select a list of Time values that are within the ranges.
Definition: timeSelector.C:88
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:759
bool processorCase() const noexcept
Return true if this is a processor case.
Definition: TimePathsI.H:29
An encapsulation of filesystem-related operations.
static mapDistributePolyMesh readProcAddressing(const faMesh &mesh, const autoPtr< faMesh > &baseMeshPtr)
Read decompose/reconstruct addressing.
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all processes in communicator.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:100
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:378
Neighbour processor patch.
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:414
void removeEmptyDirs(const fileName &path)
Remove empty directories from bottom up.
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:860
static void allGatherList(List< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses linear/tree communication. ...
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
Foam::word regionName(Foam::polyMesh::defaultRegion)
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:62
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1020
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition: topoSet.C:618
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
static void gatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
Gather data, but keep individual values separate. Uses the specified communication schedule...
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:848
void setSize(const label n)
Alias for resize()
Definition: List.H:289
static int verbose_
Output verbosity when writing.
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:671
void removeFiles(const fileName &instanceDir) const
Remove all files from mesh instance.
Definition: polyMesh.C:1330
Finite volume reconstructor for volume and surface fields.
void redistributeLagrangian(autoPtr< parLagrangianDistributor > &distributorPtr, const fvMesh &mesh, const label nOldCells, const mapDistributePolyMesh &distMap, PtrList< unmappedPassivePositionParticleCloud > &clouds)
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:584
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
const Time & time() const noexcept
Return time registry.
int dryRun() const noexcept
Return the dry-run flag.
Definition: argListI.H:109
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:397
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:118
List helper to append y unique elements onto the end of x.
Definition: ListOps.H:689
const word & system() const noexcept
Return system name.
Definition: TimePathsI.H:95
static int verbose_
Output verbosity when writing.
Reading is optional [identical to LAZY_READ].
Simple container to manage read/write, redistribute finiteArea fields.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static autoPtr< fvMesh > newMesh(const IOobject &io, const bool masterOnlyReading, const bool verbose=false)
Read mesh or create dummy mesh (0 cells, >0 patches).
Definition: fvMeshTools.C:446
Miscellaneous file handling for meshes.
static const word null
An empty word.
Definition: word.H:84
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
void removeEmptyDir(const fileName &path)
Remove empty directory.
Finite area area (element) fields.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1300
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times (as select0) otherwise return just the current ti...
Definition: timeSelector.C:265
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:997
Abstract base class for domain decomposition.
MeshObject wrapper of decompositionMethod.
static autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io, const bool decompose, const bool verbose=false)
Read mesh if available, or create empty mesh with non-proc as per proc0 mesh.
Definition: fvMeshTools.C:1183
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: argListI.H:87
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1069
static void combineReduce(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors...
bool rmDir(const fileName &directory, const bool silent=false, const bool emptyOnly=false)
Remove a directory and its contents recursively,.
Definition: POSIX.C:1433
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
void reconstructLagrangian(autoPtr< parLagrangianDistributor > &distributorPtr, const fvMesh &baseMesh, const fvMesh &mesh, const mapDistributePolyMesh &distMap, const wordRes &selectedFields)
static void writeProcAddressing(const fvMesh &procMesh, const mapDistributePolyMesh &map, const bool decompose, refPtr< fileOperation > &writeHandler)
Write addressing if decomposing (1 to many) or reconstructing (many to 1)
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:520
virtual bool write(const bool writeOnProc=true) const
Write mesh.
Definition: faMesh.C:1081
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:770
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:192
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
int debug
Static debugging option.
const fileName & rootPath() const noexcept
Return root path.
Definition: argListI.H:56
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:75
static autoPtr< faMesh > loadOrCreateMesh(const IOobject &io, const polyMesh &pMesh, const bool decompose, const bool verbose=false)
Definition: faMeshTools.C:580
static void broadcasts(const label comm, Type &arg1, Args &&... args)
Broadcast multiple items to all processes in communicator.
fileName path() const
Return the full path to the (processor local) case.
Definition: argListI.H:74
int neighbProcNo() const noexcept
Return neighbour processor number.
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))
void removeFiles(const fileName &instanceDir) const
Remove all files from mesh instance.
Definition: faMesh.C:787
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name...
Definition: instant.H:53
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:703
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: faMesh.C:747
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
static instantList findTimes(const fileName &directory, const word &constantName="constant")
Search a given directory for valid time directories.
Definition: TimePaths.C:133
static void unregisterMesh(const faMesh &mesh)
Unregister the faMesh from its associated polyMesh to prevent triggering on polyMesh changes etc...
Definition: faMeshTools.C:33
PtrList< unmappedPassivePositionParticleCloud > readLagrangian(const fvMesh &mesh, const wordList &cloudNames, const wordRes &selectedFields)
instantList times() const
Search the case for valid time directories.
Definition: TimePaths.C:142
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:667
#define WarningInFunction
Report a warning using Foam::Warning.
const fileName & globalCaseName() const noexcept
Return global case name.
Definition: TimePathsI.H:49
static List< int > & procID(const label communicator)
The list of ranks within a given communicator.
Definition: UPstream.H:1085
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: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1037
Mesh data needed to do the Finite Area discretisation.
Definition: areaFaMesh.H:47
label nPointFields
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dictionary & controlDict() const
Return read access to the controlDict dictionary.
Definition: Time.H:462
messageStream Info
Information stream (stdout output on master, null elsewhere)
void removeProcAddressing(const faMesh &mesh)
Remove procAddressing.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition: hexRef8Data.H:56
T * get() noexcept
Return pointer to managed object without nullptr checking.
Definition: autoPtr.H:216
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:89
static const fileOperation & fileHandler()
Return the current file handler. Will create the default file handler if necessary.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
UPtrList< const regIOobject > sorted() const
Return sorted list of objects.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
static autoPtr< mapDistributePolyMesh > readProcAddressing(const fvMesh &procMesh, const autoPtr< fvMesh > &baseMeshPtr)
Read procAddressing components (reconstructing)
virtual fileName objectPath(const IOobject &io, const word &typeName) const
Generate disk file name for object. Opposite of filePath.
static void forceDemandDriven(faMesh &mesh)
Force creation of everything that might vaguely be used by patches.
Definition: faMeshTools.C:45
static void setBasicGeometry(fvMesh &mesh)
Set the fvGeometryScheme to basic (to avoid parallel communication)
Definition: fvMeshTools.C:425
int verbose() const noexcept
Return the verbose flag.
Definition: argListI.H:121
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static mapDistributePolyMesh distribute(const faMesh &oldMesh, const mapDistributePolyMesh &distMap, const polyMesh &tgtPolyMesh, autoPtr< faMesh > &newMeshPtr)
Distribute mesh according to the given (volume) mesh distribution.
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
Namespace for OpenFOAM.
fileName globalPath() const
Return the full path to the global case.
Definition: argListI.H:80
static autoPtr< fileOperation > NewUncollated()
The commonly used uncollatedFileOperation.
label localSize() const
My local size.
Definition: globalIndexI.H:224
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133