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 using namespace Foam;
110 
111 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112 
113 // Use -verbose -verbose to switch on debug info. TBD.
114 int debug(::Foam::debug::debugSwitch("redistributePar", 0));
115 #define InfoOrPout (::debug ? Pout : Info.stream())
116 
117 
118 // Allocate a new file handler on valid processors only
119 // retaining the original IO ranks if possible
121 getNewHandler(const boolUList& useProc, const bool verbose = true)
122 {
123  autoPtr<fileOperation> handler
124  (
125  fileOperation::New(fileHandler(), useProc, verbose)
126  );
127 
128  if (::debug && handler)
129  {
130  Pout<< "Allocated " << handler().info()
131  << " myProcNo:" << UPstream::myProcNo(handler().comm())
132  << " ptr:" << Foam::name(handler.get()) << endl;
133  }
134 
135  return handler;
136 }
137 
138 
139 // Allocate a new file handler on valid processors only
140 // retaining the original IO ranks if possible
141 void newHandler(const boolUList& useProc, refPtr<fileOperation>& handler)
142 {
143  if (!handler)
144  {
145  handler = getNewHandler(useProc);
146  }
147 }
148 
149 
150 void createTimeDirs(const fileName& path)
151 {
152  // Get current set of local processor's time directories. Uses
153  // fileHandler
154  instantList localTimeDirs(Time::findTimes(path, "constant"));
155 
156  instantList masterTimeDirs;
157  if (Pstream::master())
158  {
159  masterTimeDirs = localTimeDirs;
160  }
161  Pstream::broadcast(masterTimeDirs);
162 
163  // Sync any cached times (e.g. masterUncollatedFileOperation::times_)
164  // since only master would have done the findTimes
165  for (const instant& t : masterTimeDirs)
166  {
167  if (!localTimeDirs.contains(t))
168  {
169  const fileName timePath(path/t.name());
170 
171  //Pout<< "Time:" << t << nl
172  // << " raw :" << timePath << nl
173  // << endl;
174  // Bypass fileHandler
175  Foam::mkDir(timePath);
176  }
177  }
178 
179  // Just to make sure remove all state and re-scan
180  fileHandler().flush();
181  (void)Time::findTimes(path, "constant");
182 }
183 
184 
185 void copyUniform
186 (
187  refPtr<fileOperation>& readHandler,
188  refPtr<fileOperation>& writeHandler,
189 
190  const bool reconstruct,
191  const bool decompose,
192 
193  const word& readTimeName,
194  const fileName& readCaseName,
195 
196  const objectRegistry& readDb,
197  const objectRegistry& writeDb
198 )
199 {
200  // 3 modes: reconstruct, decompose, redistribute
201 
202  // In reconstruct mode (separate reconstructed mesh):
203  // - read using readDb + readHandler
204  // - write using writeDb + writeHandler
205 
206  // In decompose mode (since re-using processor0 mesh):
207  // - read using readDb + readCaseName + readHandler
208  // - write using writeDb + writeHandler
209 
210  // In redistribute mode:
211  // - read using readDb + readHandler
212  // - write using writeDb + writeHandler
213 
214  fileName readPath;
215 
216  if (readHandler)
217  {
218  const label oldComm = UPstream::commWorld(readHandler().comm());
219 
220  Time& readTime = const_cast<Time&>(readDb.time());
221  bool oldProcCase = readTime.processorCase();
222  string oldCaseName;
223  if (decompose)
224  {
225  //Pout<< "***Setting caseName to " << readCaseName
226  // << " to read undecomposed uniform" << endl;
227  oldCaseName = readTime.caseName();
228  readTime.caseName() = readCaseName;
229  oldProcCase = readTime.processorCase(false);
230  }
231 
232  // Detect uniform/ at original database + time
233  readPath = readHandler().dirPath
234  (
235  false, // local directory
236  IOobject("uniform", readTimeName, readDb),
237  false // do not search in time
238  );
239 
240 
241  UPstream::commWorld(oldComm);
242 
243  if (decompose)
244  {
245  // Reset caseName on master
246  //Pout<< "***Restoring caseName to " << oldCaseName << endl;
247  readTime.caseName() = oldCaseName;
248  readTime.processorCase(oldProcCase);
249  }
250  }
252 
253  if (!readPath.empty())
254  {
255  InfoOrPout
256  << "Detected additional non-decomposed files in "
257  << readPath << endl;
258 
259  // readPath: searching is the same for all file handlers. Typical:
260  // <case>/0.1/uniform (parent dir, decompose mode)
261  // <case>/processor1/0.1/uniform (redistribute/reconstruct mode)
262  // <case>/processors2/0.1/uniform ,,
263  // writePath:
264  // uncollated : <case>/0.1/uniform (reconstruct mode). Should only
265  // be done by master
266  // uncollated : <case>/processorXXX/0.1/uniform. Should be done by all.
267  // collated : <case>/processors2/0.1/uniform. Should be done by
268  // local master only.
269 
270  const IOobject writeIO
271  (
272  "uniform",
273  writeDb.time().timeName(),
274  writeDb
275  );
276 
277  // Switch to writeHandler
278  if (writeHandler)
279  {
280  auto oldHandler = fileOperation::fileHandler(writeHandler);
281 
282  // Check: fileHandler.comm() is size 1 for uncollated
283  const label writeComm = fileHandler().comm();
284 
285  if (reconstruct)
286  {
287  const bool oldParRun = UPstream::parRun(false);
288  const label oldNumProcs(fileHandler().nProcs());
289  const fileName writePath
290  (
291  fileHandler().objectPath
292  (
293  writeIO,
294  word::null
295  )
296  );
297  fileHandler().cp(readPath, writePath);
298  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
299  UPstream::parRun(oldParRun);
300  }
301  else
302  {
303  const fileName writePath
304  (
305  fileHandler().objectPath
306  (
307  writeIO,
308  word::null
309  )
310  );
311 
312  if (::debug)
313  {
314  Pout<< " readPath :" << readPath << endl;
315  Pout<< " writePath :" << writePath << endl;
316  }
317 
318  fileHandler().broadcastCopy
319  (
320  writeComm, // send to all in writeComm
321  UPstream::master(writeComm), // to use ioranks. Check!
322  readPath,
323  writePath
324  );
325  }
326  writeHandler = fileOperation::fileHandler(oldHandler);
327  }
328  }
329 }
330 
331 
332 void printMeshData(const polyMesh& mesh)
333 {
334  // Collect all data on master
335 
336  labelListList patchNeiProcNo(Pstream::nProcs());
337  labelListList patchSize(Pstream::nProcs());
338  const labelList& pPatches = mesh.globalData().processorPatches();
339  patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
340  patchSize[Pstream::myProcNo()].setSize(pPatches.size());
341  forAll(pPatches, i)
342  {
343  const processorPolyPatch& ppp = refCast<const processorPolyPatch>
344  (
345  mesh.boundaryMesh()[pPatches[i]]
346  );
347  patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
348  patchSize[Pstream::myProcNo()][i] = ppp.size();
349  }
350  Pstream::gatherList(patchNeiProcNo);
351  Pstream::gatherList(patchSize);
352 
353 
354  // Print stats
355 
356  const globalIndex globalCells(mesh.nCells());
357  const globalIndex globalBoundaryFaces(mesh.nBoundaryFaces());
358 
359  label maxProcCells = 0;
360  label maxProcFaces = 0;
361  label totProcFaces = 0;
362  label maxProcPatches = 0;
363  label totProcPatches = 0;
364 
365  for (const int proci : Pstream::allProcs())
366  {
367  const label nLocalCells = globalCells.localSize(proci);
368  const label nBndFaces = globalBoundaryFaces.localSize(proci);
369 
370  InfoOrPout<< nl
371  << "Processor " << proci;
372 
373  if (!nLocalCells)
374  {
375  InfoOrPout<< " (empty)" << endl;
376  continue;
377  }
378  else
379  {
380  InfoOrPout<< nl
381  << " Number of cells = " << nLocalCells << endl;
382  }
383 
384  label nProcFaces = 0;
385  const labelList& nei = patchNeiProcNo[proci];
386 
387  forAll(patchNeiProcNo[proci], i)
388  {
389  InfoOrPout
390  << " Number of faces shared with processor "
391  << patchNeiProcNo[proci][i] << " = "
392  << patchSize[proci][i] << nl;
393 
394  nProcFaces += patchSize[proci][i];
395  }
396 
397  {
398  InfoOrPout
399  << " Number of processor patches = " << nei.size() << nl
400  << " Number of processor faces = " << nProcFaces << nl
401  << " Number of boundary faces = "
402  << nBndFaces-nProcFaces << endl;
403  }
404 
405  maxProcCells = max(maxProcCells, nLocalCells);
406  totProcFaces += nProcFaces;
407  totProcPatches += nei.size();
408  maxProcFaces = max(maxProcFaces, nProcFaces);
409  maxProcPatches = max(maxProcPatches, nei.size());
410  }
411 
412  // Summary stats
413 
414  InfoOrPout
415  << nl
416  << "Number of processor faces = " << (totProcFaces/2) << nl
417  << "Max number of cells = " << maxProcCells;
418 
419  if (maxProcCells != globalCells.totalSize())
420  {
421  scalar avgValue = scalar(globalCells.totalSize())/Pstream::nProcs();
422 
423  InfoOrPout
424  << " (" << 100.0*(maxProcCells-avgValue)/avgValue
425  << "% above average " << avgValue << ')';
426  }
427  InfoOrPout<< nl;
428 
429  InfoOrPout<< "Max number of processor patches = " << maxProcPatches;
430  if (totProcPatches)
431  {
432  scalar avgValue = scalar(totProcPatches)/Pstream::nProcs();
433 
434  InfoOrPout
435  << " (" << 100.0*(maxProcPatches-avgValue)/avgValue
436  << "% above average " << avgValue << ')';
437  }
438  InfoOrPout<< nl;
439 
440  InfoOrPout<< "Max number of faces between processors = " << maxProcFaces;
441  if (totProcFaces)
442  {
443  scalar avgValue = scalar(totProcFaces)/Pstream::nProcs();
444 
445  InfoOrPout
446  << " (" << 100.0*(maxProcFaces-avgValue)/avgValue
447  << "% above average " << avgValue << ')';
448  }
449  InfoOrPout<< nl << endl;
450 }
451 
452 
453 // Debugging: write volScalarField with decomposition for post processing.
454 void writeDecomposition
455 (
456  const word& name,
457  const fvMesh& mesh,
458  const labelUList& decomp
459 )
460 {
461  // Write the decomposition as labelList for use with 'manual'
462  // decomposition method.
464  (
465  IOobject
466  (
467  "cellDecomposition",
468  mesh.facesInstance(), // mesh read from facesInstance
469  mesh,
473  ),
474  decomp
475  ).write();
476 
477  InfoOrPout
478  << "Writing wanted cell distribution to volScalarField " << name
479  << " for postprocessing purposes." << nl << endl;
480 
481  volScalarField procCells
482  (
483  IOobject
484  (
485  name,
486  mesh.time().timeName(),
487  mesh,
491  ),
492  mesh,
493  dimensionedScalar("", dimless, -1),
495  );
496 
497  forAll(procCells, celli)
498  {
499  procCells[celli] = decomp[celli];
500  }
501 
502  procCells.correctBoundaryConditions();
503  procCells.write();
504 }
505 
506 
507 void determineDecomposition
508 (
509  refPtr<fileOperation>& readHandler,
510  const Time& baseRunTime,
511  const fileName& decompDictFile, // optional location for decomposeParDict
512  const bool decompose, // decompose, i.e. read from undecomposed case
513  const fileName& proc0CaseName,
514  const fvMesh& mesh,
515  const bool writeCellDist,
516 
517  label& nDestProcs,
518  labelList& decomp
519 )
520 {
521  // Switch to readHandler since decomposition method might do IO
522  // (e.g. read decomposeParDict)
523  auto oldHandler = fileOperation::fileHandler(readHandler);
524 
525  // Read decomposeParDict (on all processors)
527  (
528  mesh,
529  decompDictFile
530  );
531 
532  decompositionMethod& decomposer = method.decomposer();
533 
534  if (!decomposer.parallelAware())
535  {
537  << "You have selected decomposition method \""
538  << decomposer.type() << "\n"
539  << " which does not synchronise decomposition across"
540  " processor patches.\n"
541  " You might want to select a decomposition method"
542  " that is aware of this. Continuing...." << endl;
543  }
544 
545  Time& tm = const_cast<Time&>(mesh.time());
546 
547  const bool oldProcCase = tm.processorCase();
548  if (decompose)
549  {
550  InfoOrPout
551  << "Setting caseName to " << baseRunTime.caseName()
552  << " to read decomposeParDict" << endl;
553  tm.caseName() = baseRunTime.caseName();
554  tm.processorCase(false);
555  }
556 
557  scalarField cellWeights;
558  if (method.found("weightField"))
559  {
560  word weightName = method.get<word>("weightField");
561 
562  volScalarField weights
563  (
564  IOobject
565  (
566  weightName,
567  tm.timeName(),
568  mesh,
572  ),
573  mesh
574  );
575  cellWeights = weights.internalField();
576  }
577 
578  nDestProcs = decomposer.nDomains();
579  decomp = decomposer.decompose(mesh, cellWeights);
580 
581  readHandler = fileOperation::fileHandler(oldHandler);
582 
583 
584  if (decompose)
585  {
586  InfoOrPout
587  << "Restoring caseName" << endl;
588  tm.caseName() = proc0CaseName;
589  tm.processorCase(oldProcCase);
590  }
591 
592  // Dump decomposition to volScalarField
593  if (writeCellDist)
594  {
595  const label oldNumProcs =
596  const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs);
597 
598  // Note: on master make sure to write to processor0
599  if (decompose)
600  {
601  if (UPstream::master())
602  {
603  const bool oldParRun = UPstream::parRun(false);
604 
605  InfoOrPout
606  << "Setting caseName to " << baseRunTime.caseName()
607  << " to write undecomposed cellDist" << endl;
608 
609  tm.caseName() = baseRunTime.caseName();
610  tm.processorCase(false);
611  writeDecomposition("cellDist", mesh, decomp);
612  InfoOrPout<< "Restoring caseName" << endl;
613  tm.caseName() = proc0CaseName;
614  tm.processorCase(oldProcCase);
615 
616  UPstream::parRun(oldParRun);
617  }
618  }
619  else
620  {
621  writeDecomposition("cellDist", mesh, decomp);
622  }
623 
624  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
625  }
626 }
627 
628 
629 // Variant of GeometricField::correctBoundaryConditions that only
630 // evaluates selected patch fields
631 template<class CoupledPatchType, class GeoField>
632 void correctCoupledBoundaryConditions(fvMesh& mesh)
633 {
634  for (GeoField& fld : mesh.sorted<GeoField>())
635  {
636  fld.boundaryFieldRef().template evaluateCoupled<CoupledPatchType>();
637  }
638 }
639 
640 
641 // Inplace redistribute mesh and any fields
642 autoPtr<mapDistributePolyMesh> redistributeAndWrite
643 (
644  refPtr<fileOperation>& readHandler,
645  refPtr<fileOperation>& writeHandler,
646  const Time& baseRunTime,
647  const fileName& proc0CaseName,
648 
649  // Controls
650  const bool doReadFields,
651  const bool decompose, // decompose, i.e. read from undecomposed case
652  const bool reconstruct,
653  const bool overwrite,
654 
655  // Decomposition information
656  const label nDestProcs,
657  const labelList& decomp,
658 
659  // Mesh information
660  const boolList& volMeshOnProc,
661  const fileName& volMeshInstance,
662  fvMesh& mesh
663 )
664 {
665  Time& runTime = const_cast<Time&>(mesh.time());
666  const bool oldProcCase = runTime.processorCase();
667 
669  //Pout<< "Before distribution:" << endl;
670  //printMeshData(mesh);
671 
672 
673  // Storage of fields
674 
675  PtrList<volScalarField> volScalarFields;
676  PtrList<volVectorField> volVectorFields;
677  PtrList<volSphericalTensorField> volSphereTensorFields;
678  PtrList<volSymmTensorField> volSymmTensorFields;
679  PtrList<volTensorField> volTensorFields;
680 
681  PtrList<surfaceScalarField> surfScalarFields;
682  PtrList<surfaceVectorField> surfVectorFields;
683  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
684  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
685  PtrList<surfaceTensorField> surfTensorFields;
686 
692 
693  PtrList<pointScalarField> pointScalarFields;
694  PtrList<pointVectorField> pointVectorFields;
695  PtrList<pointTensorField> pointTensorFields;
696  PtrList<pointSphericalTensorField> pointSphTensorFields;
697  PtrList<pointSymmTensorField> pointSymmTensorFields;
698 
699  // Self-contained pointMesh for reading pointFields
700  const pointMesh oldPointMesh(mesh);
701 
702  // Track how many (if any) pointFields are read/mapped
703  label nPointFields = 0;
704 
705  refPtr<fileOperation> noWriteHandler;
706 
707  parPointFieldDistributor pointDistributor
708  (
709  oldPointMesh, // source mesh
710  false, // savePoints=false (ie, delay until later)
711  //false // Do not write
712  noWriteHandler // Do not write
713  );
714 
715 
716  if (doReadFields)
717  {
718  // Create 0 sized mesh to do all the generation of zero sized
719  // fields on processors that have zero sized meshes. Note that this is
720  // only necessary on master but since polyMesh construction with
721  // Pstream::parRun does parallel comms we have to do it on all
722  // processors
723  autoPtr<fvMeshSubset> subsetterPtr;
724 
725  // Missing a volume mesh somewhere?
726  if (volMeshOnProc.found(false))
727  {
728  // A zero-sized mesh with boundaries.
729  // This is used to create zero-sized fields.
730  subsetterPtr.reset(new fvMeshSubset(mesh, zero{}));
731  subsetterPtr().subMesh().init(true);
732  subsetterPtr().subMesh().globalData();
733  subsetterPtr().subMesh().tetBasePtIs();
734  subsetterPtr().subMesh().geometricD();
735  }
736 
737 
738  // Get original objects (before incrementing time!)
739  if (decompose)
740  {
741  InfoOrPout
742  << "Setting caseName to " << baseRunTime.caseName()
743  << " to read IOobjects" << endl;
744  runTime.caseName() = baseRunTime.caseName();
745  runTime.processorCase(false);
746  }
747 
748  //IOobjectList objects(mesh, runTime.timeName());
749  // Swap to reading fileHandler and read IOobjects
750  IOobjectList objects;
751  if (readHandler)
752  {
753  auto oldHandler = fileOperation::fileHandler(readHandler);
754  const label oldComm = UPstream::commWorld(fileHandler().comm());
755 
756  objects = IOobjectList(mesh, runTime.timeName());
757  readHandler = fileOperation::fileHandler(oldHandler);
758  UPstream::commWorld(oldComm);
759  }
760 
761  if (decompose)
762  {
763  InfoOrPout
764  << "Restoring caseName" << endl;
765  runTime.caseName() = proc0CaseName;
766  runTime.processorCase(oldProcCase);
767  }
768 
769  InfoOrPout
770  << "From time " << runTime.timeName()
771  << " mesh:" << mesh.objectRegistry::objectRelPath()
772  << " have objects:" << objects.names() << endl;
773 
774  // We don't want to map the decomposition (mapping already tested when
775  // mapping the cell centre field)
776  (void)objects.erase("cellDist");
777 
778 
779  if (decompose)
780  {
781  runTime.caseName() = baseRunTime.caseName();
782  runTime.processorCase(false);
783  }
784 
785  // Field reading
786 
787  #undef doFieldReading
788  #define doFieldReading(Storage) \
789  { \
790  fieldsDistributor::readFields \
791  ( \
792  volMeshOnProc, readHandler, mesh, subsetterPtr, \
793  objects, Storage \
794  ); \
795  }
796 
797  // volField
798  doFieldReading(volScalarFields);
799  doFieldReading(volVectorFields);
800  doFieldReading(volSphereTensorFields);
801  doFieldReading(volSymmTensorFields);
802  doFieldReading(volTensorFields);
803 
804  // surfaceField
805  doFieldReading(surfScalarFields);
806  doFieldReading(surfVectorFields);
807  doFieldReading(surfSphereTensorFields);
808  doFieldReading(surfSymmTensorFields);
809  doFieldReading(surfTensorFields);
810 
811  // Dimensioned internal fields
812  doFieldReading(dimScalarFields);
813  doFieldReading(dimVectorFields);
814  doFieldReading(dimSphereTensorFields);
815  doFieldReading(dimSymmTensorFields);
816  doFieldReading(dimTensorFields);
817 
818  // pointFields
819  nPointFields = 0;
820 
821  #undef doFieldReading
822  #define doFieldReading(Storage) \
823  { \
824  fieldsDistributor::readFields \
825  ( \
826  volMeshOnProc, readHandler, oldPointMesh, \
827  subsetterPtr, objects, Storage, \
828  true /* (deregister field) */ \
829  ); \
830  nPointFields += Storage.size(); \
831  }
832 
833  doFieldReading(pointScalarFields);
834  doFieldReading(pointVectorFields);
835  doFieldReading(pointSphTensorFields);
836  doFieldReading(pointSymmTensorFields);
837  doFieldReading(pointTensorFields);
838  #undef doFieldReading
839 
840 
841  // Done reading
842 
843  if (decompose)
844  {
845  runTime.caseName() = proc0CaseName;
846  runTime.processorCase(oldProcCase);
847  }
848  }
849 
850  // Save pointMesh information before any topology changes occur!
851  if (nPointFields)
852  {
853  pointDistributor.saveMeshPoints();
854  }
855 
856 
857  // Read refinement data
858  autoPtr<hexRef8Data> refDataPtr;
859  {
860  // Read refinement data
861  if (decompose)
862  {
863  runTime.caseName() = baseRunTime.caseName();
864  runTime.processorCase(false);
865  }
866  IOobject io
867  (
868  "dummy",
869  volMeshInstance, //mesh.facesInstance(),
871  mesh,
875  );
876 
877  if (readHandler)
878  {
879  auto oldHandler = fileOperation::fileHandler(readHandler);
880  const label oldComm = UPstream::commWorld(fileHandler().comm());
881 
882  // Read
883  refDataPtr.reset(new hexRef8Data(io));
884 
885  UPstream::commWorld(oldComm);
886  readHandler = fileOperation::fileHandler(oldHandler);
887  }
888  else
889  {
891  refDataPtr.reset(new hexRef8Data(io));
892  }
893 
894  if (decompose)
895  {
896  runTime.caseName() = proc0CaseName;
897  runTime.processorCase(oldProcCase);
898  }
899 
900  // Make sure all processors have valid data (since only some will
901  // read)
902  refDataPtr().sync(io);
903 
904  // Now we've read refinement data we can remove it
906  }
907 
908 
909  // Mesh distribution engine
910  fvMeshDistribute distributor(mesh);
911 
912  // Do all the distribution of mesh and fields
913  autoPtr<mapDistributePolyMesh> distMap = distributor.distribute(decomp);
914 
915  // Print some statistics
916  InfoOrPout<< "After distribution:" << endl;
917  printMeshData(mesh);
918 
919  // Get other side of processor boundaries
920  do
921  {
922  #undef doCorrectCoupled
923  #define doCorrectCoupled(FieldType) \
924  correctCoupledBoundaryConditions<processorFvPatch, FieldType>(mesh);
925 
926  doCorrectCoupled(volScalarField);
927  doCorrectCoupled(volVectorField);
928  doCorrectCoupled(volSphericalTensorField);
929  doCorrectCoupled(volSymmTensorField);
930  doCorrectCoupled(volTensorField);
931  #undef doCorrectCoupled
932  }
933  while (false);
934 
935  // No update surface fields
936 
937 
938  // Map pointFields
939  if (nPointFields)
940  {
941  // Construct new pointMesh from distributed mesh
942  const pointMesh& newPointMesh = pointMesh::New(mesh);
943 
944  pointDistributor.resetTarget(newPointMesh, distMap());
945 
946  pointDistributor.distributeAndStore(pointScalarFields);
947  pointDistributor.distributeAndStore(pointVectorFields);
948  pointDistributor.distributeAndStore(pointSphTensorFields);
949  pointDistributor.distributeAndStore(pointSymmTensorFields);
950  pointDistributor.distributeAndStore(pointTensorFields);
951  }
952 
953 
954  // Set the minimum write precision
956 
957 
958  if (!overwrite)
959  {
960  ++runTime;
962  }
963  else
964  {
965  mesh.setInstance(volMeshInstance);
966  }
967 
968 
969  // Register mapDistributePolyMesh for automatic writing...
970  IOmapDistributePolyMeshRef distMapRef
971  (
972  IOobject
973  (
974  "procAddressing",
977  mesh.thisDb(),
980  ),
981  distMap()
982  );
983 
984 
985  if (reconstruct)
986  {
987  auto oldHandler = fileOperation::fileHandler(writeHandler);
988 
989  if (UPstream::master())
990  {
991  InfoOrPout
992  << "Setting caseName to " << baseRunTime.caseName()
993  << " to write reconstructed mesh (and fields)." << endl;
994  runTime.caseName() = baseRunTime.caseName();
995  const bool oldProcCase(runTime.processorCase(false));
996  const label oldNumProcs
997  (
998  const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs)
999  );
1000  const bool oldParRun = UPstream::parRun(false);
1001 
1002  mesh.write();
1004 
1005  // Now we've written all. Reset caseName on master
1006  InfoOrPout<< "Restoring caseName" << endl;
1007  UPstream::parRun(oldParRun);
1008  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
1009  runTime.caseName() = proc0CaseName;
1010  runTime.processorCase(oldProcCase);
1011  }
1012 
1013  writeHandler = fileOperation::fileHandler(oldHandler);
1014  }
1015  else
1016  {
1017  auto oldHandler = fileOperation::fileHandler(writeHandler);
1018 
1019  const label oldNumProcs
1020  (
1021  const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs)
1022  );
1023  mesh.write();
1024  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
1025 
1026  writeHandler = fileOperation::fileHandler(oldHandler);
1027 
1029  }
1030  InfoOrPout
1031  << "Written redistributed mesh to "
1032  << mesh.facesInstance() << nl << endl;
1033 
1034 
1035  if (decompose)
1036  {
1037  // Decompose (1 -> N)
1038  // so {boundary,cell,face,point}ProcAddressing have meaning
1040  (
1041  mesh,
1042  distMap(),
1043  decompose,
1044  mesh.facesInstance(), //oldFacesInstance,
1045  writeHandler // to write *ProcAddressing
1046  );
1047  }
1048  else if (reconstruct)
1049  {
1050  // Reconstruct (N -> 1)
1051  // so {boundary,cell,face,point}ProcAddressing have meaning. Make sure
1052  // to write these to meshes containing the source meshes (i.e. using
1053  // the read handler)
1055  (
1056  mesh,
1057  distMap(),
1058  decompose,
1059  volMeshInstance, //oldFacesInstance,
1060  readHandler //writeHandler
1061  );
1062  }
1063  else
1064  {
1065  // Redistribute (N -> M)
1066  // {boundary,cell,face,point}ProcAddressing would be incorrect
1067  // - can either remove or redistribute previous
1069  }
1070 
1071 
1072  // Refinement data
1073  if (refDataPtr)
1074  {
1075  auto& refData = refDataPtr();
1076 
1077  // Set instance
1078  IOobject io
1079  (
1080  "dummy",
1081  mesh.facesInstance(),
1083  mesh,
1087  );
1088  refData.sync(io);
1089 
1090 
1091  // Distribute
1092  refData.distribute(distMap());
1093 
1094 
1095  auto oldHandler = fileOperation::fileHandler(writeHandler);
1096 
1097  if (reconstruct)
1098  {
1099  if (UPstream::master())
1100  {
1101  const bool oldParRun = UPstream::parRun(false);
1102  const label oldNumProcs
1103  (
1104  const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs)
1105  );
1106 
1107  InfoOrPout
1108  << "Setting caseName to " << baseRunTime.caseName()
1109  << " to write reconstructed refinement data." << endl;
1110  runTime.caseName() = baseRunTime.caseName();
1111  const bool oldProcCase(runTime.processorCase(false));
1112 
1113  refData.write();
1114 
1115  // Now we've written all. Reset caseName on master
1116  InfoOrPout<< "Restoring caseName" << endl;
1117  runTime.caseName() = proc0CaseName;
1118  runTime.processorCase(oldProcCase);
1119 
1120  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
1121  UPstream::parRun(oldParRun);
1122  }
1123  }
1124  else
1125  {
1126  const label oldNumProcs
1127  (
1128  const_cast<fileOperation&>(fileHandler()).nProcs(nDestProcs)
1129  );
1130  refData.write();
1131  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
1132  }
1133 
1134  writeHandler = fileOperation::fileHandler(oldHandler);
1135  }
1136 
1138  //{
1139  // // Read sets
1140  // if (decompose)
1141  // {
1142  // runTime.caseName() = baseRunTime.caseName();
1143  // runTime.processorCase(false);
1144  // }
1145  // IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1146  //
1147  // PtrList<cellSet> cellSets;
1148  // ReadFields(objects, cellSets);
1149  //
1150  // if (decompose)
1151  // {
1152  // runTime.caseName() = proc0CaseName;
1153  // runTime.processorCase(oldProcCase);
1154  // }
1155  //
1156  // forAll(cellSets, i)
1157  // {
1158  // cellSets[i].distribute(distMap());
1159  // }
1160  //
1161  // if (reconstruct)
1162  // {
1163  // if (Pstream::master())
1164  // {
1165  // Pout<< "Setting caseName to " << baseRunTime.caseName()
1166  // << " to write reconstructed refinement data." << endl;
1167  // runTime.caseName() = baseRunTime.caseName();
1168  // const bool oldProcCase(runTime.processorCase(false));
1169  //
1170  // forAll(cellSets, i)
1171  // {
1172  // cellSets[i].distribute(distMap());
1173  // }
1174  //
1175  // // Now we've written all. Reset caseName on master
1176  // Pout<< "Restoring caseName" << endl;
1177  // runTime.caseName() = proc0CaseName;
1178  // runTime.processorCase(oldProcCase);
1179  // }
1180  // }
1181  // else
1182  // {
1183  // forAll(cellSets, i)
1184  // {
1185  // cellSets[i].distribute(distMap());
1186  // }
1187  // }
1188  //}
1189 
1190 
1191  return distMap;
1192 }
1193 
1194 
1195 /*---------------------------------------------------------------------------*\
1196  main
1197 \*---------------------------------------------------------------------------*/
1198 
1199 int main(int argc, char *argv[])
1200 {
1202  (
1203  "Redistribute decomposed mesh and fields according"
1204  " to the decomposeParDict settings.\n"
1205  "Optionally run in decompose/reconstruct mode"
1206  );
1207 
1208  argList::noFunctionObjects(); // Never use function objects
1209 
1210  // enable -constant ... if someone really wants it
1211  // enable -zeroTime to prevent accidentally trashing the initial fields
1212  timeSelector::addOptions(true, true);
1213 
1214  #include "addAllRegionOptions.H"
1215 
1216  #include "addOverwriteOption.H"
1217  argList::addBoolOption("decompose", "Decompose case");
1218  argList::addBoolOption("reconstruct", "Reconstruct case");
1220  (
1221  "Additional verbosity. (Can be used multiple times for debugging)"
1222  );
1224  (
1225  "Test without writing the decomposition. "
1226  "Changes -cellDist to only write volScalarField."
1227  );
1228  argList::addVerboseOption("Additional verbosity");
1230  (
1231  "cellDist",
1232  "Write cell distribution as a labelList - for use with 'manual' "
1233  "decomposition method or as a volScalarField for post-processing."
1234  );
1236  (
1237  "newTimes",
1238  "Only reconstruct new times (i.e. that do not exist already)"
1239  );
1241  (
1242  "Additional verbosity. (Can be used multiple times)"
1243  );
1245  (
1246  "no-finite-area",
1247  "Suppress finiteArea mesh/field handling",
1248  true // Advanced option
1249  );
1250 
1251 
1252  //- Disable caching of times/processor dirs etc. Cause massive parallel
1253  // problems when e.g decomposing.
1255 
1256 
1257  // Handle arguments
1258  // ~~~~~~~~~~~~~~~~
1259  // (replacement for setRootCase that does not abort)
1260 
1261  argList args(argc, argv);
1262 
1263 
1264  // As much as possible avoid synchronised operation. To be looked at more
1265  // closely for the three scenarios:
1266  // - decompose - reads on master (and from parent directory) and sends
1267  // dictionary to slaves
1268  // - distribute - reads on potentially a different number of processors
1269  // than it writes to
1270  // - reconstruct - reads parallel, write on master only and to parent
1271  // directory
1272 
1273  // Running data-distributed?
1274  // (processors cannot see all other processors' files)
1275  const bool hasDistributedFiles(fileHandler().distributed());
1276 
1277 
1278  // Set up loose processorsXXX directory matching (for collated) so e.g.
1279  // when checking for -latestTime we don't miss anything. Once we know
1280  // the time, actual number of processors etc we switch back to strict
1281  // matching.
1283 
1284  // Need this line since we don't include "setRootCase.H"
1285  #include "foamDlOpenLibs.H"
1286 
1287  const bool reconstruct = args.found("reconstruct");
1288  const bool writeCellDist = args.found("cellDist");
1289  const bool dryrun = args.dryRun();
1290  const bool newTimes = args.found("newTimes");
1291  const int optVerbose = args.verbose();
1292 
1293  const bool doFiniteArea = !args.found("no-finite-area");
1294  bool decompose = args.found("decompose");
1295  bool overwrite = args.found("overwrite");
1296 
1297  // Field restrictions...
1298  const wordRes selectedFields;
1299  const wordRes selectedLagrangianFields;
1300 
1301 
1302  if (!UPstream::parRun())
1303  {
1305  << ": This utility can only be run parallel"
1306  << exit(FatalError);
1307  }
1308 
1309  if (decompose && reconstruct)
1310  {
1312  << "Cannot specify both -decompose and -reconstruct"
1313  << exit(FatalError);
1314  }
1315 
1316  if (optVerbose)
1317  {
1318  // Report on output
1321 
1322  if (optVerbose > 1)
1323  {
1324  Info<< "Additional debugging enabled" << nl << endl;
1325  ::debug = 1;
1326  }
1327  }
1328 
1329  // Disable NaN setting and floating point error trapping. This is to avoid
1330  // any issues inside the field redistribution inside fvMeshDistribute
1331  // which temporarily moves processor faces into existing patches. These
1332  // will now not have correct values. After all bits have been assembled
1333  // the processor fields will be restored but until then there might
1334  // be uninitialised values which might upset some patch field constructors.
1335  // Workaround by disabling floating point error trapping. TBD: have
1336  // actual field redistribution instead of subsetting inside
1337  // fvMeshDistribute.
1338  Foam::sigFpe::unset(true);
1339 
1340 
1341  // File handlers (read/write)
1342  // ==========================
1343 
1344  // Read handler on processors with a volMesh
1345  refPtr<fileOperation> volMeshReadHandler;
1346 
1347  // Read handler on processors with an areaMesh
1348  refPtr<fileOperation> areaMeshReadHandler;
1349 
1350  // Handler for master-only operation (read/writing from/to undecomposed)
1351  // - only the 'real' master, not io-rank masters
1352  refPtr<fileOperation> masterOnlyHandler;
1353 
1355  {
1356  const bool oldParRun = UPstream::parRun(false);
1357 
1358  masterOnlyHandler = fileOperation::NewUncollated();
1359 
1360  UPstream::parRun(oldParRun);
1361  }
1362 
1363  if (decompose)
1364  {
1365  InfoOrPout<< "Decomposing case (like decomposePar)"
1366  << nl << endl;
1367  }
1368  else if (reconstruct)
1369  {
1370  InfoOrPout<< "Reconstructing case (like reconstructParMesh)"
1371  << nl << endl;
1372  }
1373 
1374  if (decompose || reconstruct)
1375  {
1376  // The UPstream::nProcs is either the source or destination procs
1378  InfoOrPout<< "Switching to exact matching for "
1380  << " processor directories"
1381  << nl << endl;
1382  }
1383  else
1384  {
1385  // Redistribute mode. Accept any processorsXXX naming since we don't
1386  // know yet what the source/target number of processors is
1388  InfoOrPout<< "Switching to matching any "
1389  << fileOperation::processorsBaseDir << " directory" << nl << endl;
1390  }
1391 
1392 
1393  if ((decompose || reconstruct) && !overwrite)
1394  {
1395  overwrite = true;
1396 
1397  Warning
1398  << nl << " "
1399  << "Added implicit -overwrite for decompose/reconstruct modes"
1400  << nl << endl;
1401  }
1402 
1403  if
1404  (
1405  fileHandler().ioRanks().contains(UPstream::myProcNo())
1406  && !Foam::isDir(args.rootPath())
1407  )
1408  {
1409  //FatalErrorInFunction
1411  << ": cannot open root directory " << args.rootPath()
1412  << endl;
1413  //<< exit(FatalError);
1414  }
1415 
1416 
1417  if (hasDistributedFiles)
1418  {
1419  InfoOrPout<< "Detected multiple roots i.e. non-nfs running"
1420  << nl << endl;
1421 
1422  fileHandler().mkDir(args.globalPath());
1423  }
1424 
1425 
1426  // Check if we have processor directories. Ideally would like to
1427  // use fileHandler().dirPath here but we don't have runTime yet and
1428  // want to delay constructing runTime until we've synced all time
1429  // directories...
1430  const fileName procDir(fileHandler().filePath(args.path()));
1431  if (Foam::isDir(procDir))
1432  {
1433  if (decompose)
1434  {
1435  InfoOrPout<< "Removing existing processor directory:"
1436  << args.relativePath(procDir) << endl;
1437  fileHandler().rmDir(procDir);
1438  }
1439  }
1440  else
1441  {
1442  // Directory does not exist. If this happens on master -> decompose mode
1443  if (UPstream::master() && !reconstruct && !decompose)
1444  {
1445  decompose = true;
1446  InfoOrPout
1447  << "No processor directories; switching on decompose mode"
1448  << nl << endl;
1449  }
1450  }
1451  // If master changed to decompose mode make sure all nodes know about it
1452  Pstream::broadcast(decompose);
1453  if (decompose)
1454  {
1455  // The UPstream::nProcs is either the source or destination procs
1457  InfoOrPout
1458  << "Switching to exact matching for "
1460  << " processor directories"
1461  << nl << endl;
1462  }
1463 
1464 
1465 
1466  // If running distributed we have problem of new processors not finding
1467  // a system/controlDict. However if we switch on the master-only reading
1468  // the problem becomes that the time directories are differing sizes and
1469  // e.g. latestTime will pick up a different time (which causes createTime.H
1470  // to abort). So for now make sure to have master times on all
1471  // processors
1472  if (!reconstruct)
1473  {
1474  InfoOrPout<< "Creating time directories on all processors"
1475  << nl << endl;
1476  createTimeDirs(args.path());
1477  }
1478 
1479  // Construct time
1480  // ~~~~~~~~~~~~~~
1481 
1482  // Replace #include "createTime.H" with our own version
1483  // that has MUST_READ instead of READ_MODIFIED
1484 
1485  Info<< "Create time\n" << endl;
1486  Time runTime
1487  (
1489  args,
1490  false, // no enableFunctionObjects
1491  true, // permit enableLibs
1492  IOobjectOption::MUST_READ // Without watching
1493  );
1494 
1495 
1496  refPtr<fileOperation> writeHandler;
1497 
1498 
1499  runTime.functionObjects().off(); // Extra safety?
1500  // Make sure that no runTime checking is done since fileOperation::addWatch
1501  // etc. does broadcast over world, even if constructed only on a subset
1502  // of procs
1503  runTime.runTimeModifiable(false);
1504 
1505  // Save local processor0 casename
1506  const fileName proc0CaseName = runTime.caseName();
1507  const bool oldProcCase = runTime.processorCase();
1508 
1509 
1510  // Construct undecomposed Time
1511  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1512  // This will read the same controlDict but might have a different
1513  // set of times so enforce same times
1514 
1515  if (hasDistributedFiles)
1516  {
1517  InfoOrPout<< "Creating time directories for undecomposed Time"
1518  << " on all processors" << nl << endl;
1519  createTimeDirs(args.globalPath());
1520  }
1521 
1522 
1523  InfoOrPout<< "Create undecomposed database" << nl << endl;
1524  Time baseRunTime
1525  (
1526  runTime.controlDict(),
1527  runTime.rootPath(),
1529  runTime.system(),
1530  runTime.constant(),
1531  false, // No function objects
1532  false, // No extra controlDict libs
1533  IOobjectOption::MUST_READ // Without watching
1534  );
1535  // Make sure that no runTime checking is done since fileOperation::addWatch
1536  // etc. does broadcast over world, even if constructed only on a subset
1537  // of procs
1538  baseRunTime.runTimeModifiable(false);
1539 
1540 
1541  wordHashSet masterTimeDirSet;
1542  if (newTimes)
1543  {
1544  instantList baseTimeDirs(baseRunTime.times());
1545  for (const instant& t : baseTimeDirs)
1546  {
1547  masterTimeDirSet.insert(t.name());
1548  }
1549  }
1550 
1551 
1552  // Allow override of decomposeParDict location
1553  const fileName decompDictFile =
1554  args.getOrDefault<fileName>("decomposeParDict", "");
1555 
1556  // Get region names
1557  #include "getAllRegionOptions.H"
1558 
1560  {
1561  InfoOrPout<< "Using region: " << regionNames[0] << nl << endl;
1562  }
1563 
1564 
1565  // Demand driven lagrangian mapper
1566  autoPtr<parLagrangianDistributor> lagrangianDistributorPtr;
1567 
1568  if (reconstruct)
1569  {
1570  // use the times list from the master processor
1571  // and select a subset based on the command-line options
1573  Pstream::broadcast(timeDirs);
1574 
1575  if (timeDirs.empty())
1576  {
1578  << "No times selected"
1579  << exit(FatalError);
1580  }
1581 
1582 
1583  // Pass1 : reconstruct mesh and addressing
1584  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1585 
1586 
1587  InfoOrPout
1588  << nl
1589  << "Reconstructing mesh and addressing" << nl << endl;
1590 
1591  forAll(regionNames, regioni)
1592  {
1593  const word& regionName = regionNames[regioni];
1594 
1595  const fileName volMeshSubDir
1596  (
1598  );
1599  const fileName areaMeshSubDir
1600  (
1602  );
1603 
1604  InfoOrPout
1605  << nl
1606  << "Reconstructing mesh:"
1608 
1609  bool areaMeshDetected = false;
1610 
1611  // Loop over all times
1612  forAll(timeDirs, timeI)
1613  {
1614  // Set time for global database
1615  runTime.setTime(timeDirs[timeI], timeI);
1616  baseRunTime.setTime(timeDirs[timeI], timeI);
1617 
1618  InfoOrPout<< "Time = " << runTime.timeName() << endl << endl;
1619 
1620  // Where meshes are
1621  fileName volMeshInstance;
1622  fileName areaMeshInstance;
1623 
1624  volMeshInstance = runTime.findInstance
1625  (
1626  volMeshSubDir,
1627  "faces",
1629  );
1630 
1631  if (doFiniteArea)
1632  {
1633  areaMeshInstance = runTime.findInstance
1634  (
1635  areaMeshSubDir,
1636  "faceLabels",
1638  );
1639  }
1640 
1642  (
1644  volMeshInstance,
1645  areaMeshInstance
1646  );
1647 
1648 
1649  // Check processors have meshes
1650  // - check for 'faces' file (polyMesh)
1651  // - check for 'faceLabels' file (faMesh)
1652  boolList volMeshOnProc;
1653  boolList areaMeshOnProc(UPstream::nProcs(), false);
1654 
1655  volMeshOnProc = haveMeshFile
1656  (
1657  runTime,
1658  volMeshInstance/volMeshSubDir,
1659  "faces"
1660  );
1661 
1662  // Create handler for reading
1663  newHandler(volMeshOnProc, volMeshReadHandler);
1664 
1665  if (doFiniteArea)
1666  {
1667  areaMeshOnProc = haveMeshFile
1668  (
1669  runTime,
1670  areaMeshInstance/areaMeshSubDir,
1671  "faceLabels"
1672  );
1673  areaMeshDetected = areaMeshOnProc.found(true);
1674 
1675  if (areaMeshOnProc == volMeshOnProc)
1676  {
1677  if (volMeshReadHandler)
1678  {
1679  // Use same reader for faMesh as for fvMesh
1680  areaMeshReadHandler.ref(volMeshReadHandler.ref());
1681  }
1682  }
1683  else
1684  {
1685  newHandler(areaMeshOnProc, areaMeshReadHandler);
1686  }
1687  }
1688 
1689 
1690  // Addressing back to reconstructed mesh as xxxProcAddressing.
1691  // - all processors have consistent faceProcAddressing
1692  // - processors without a mesh don't need faceProcAddressing
1693 
1694 
1695  // Note: filePath searches up on processors that don't have
1696  // processor if instance = constant so explicitly check
1697  // found filename.
1698  bool haveVolAddressing = false;
1699  if (volMeshOnProc[Pstream::myProcNo()])
1700  {
1701  // Read faces (just to know their size)
1702  faceCompactIOList faces
1703  (
1704  IOobject
1705  (
1706  "faces",
1707  volMeshInstance,
1708  volMeshSubDir,
1709  runTime,
1711  )
1712  );
1713 
1714  // Check faceProcAddressing
1715  labelIOList faceProcAddressing
1716  (
1717  IOobject
1718  (
1719  "faceProcAddressing",
1720  volMeshInstance,
1721  volMeshSubDir,
1722  runTime,
1724  )
1725  );
1726 
1727  haveVolAddressing =
1728  (
1729  faceProcAddressing.headerOk()
1730  && faceProcAddressing.size() == faces.size()
1731  );
1732  }
1733  else
1734  {
1735  // Have no mesh. Don't need addressing
1736  haveVolAddressing = true;
1737  }
1738 
1739  bool haveAreaAddressing = false;
1740  if (areaMeshOnProc[Pstream::myProcNo()])
1741  {
1742  // Read faces (just to know their size)
1744  (
1745  IOobject
1746  (
1747  "faceLabels",
1748  areaMeshInstance,
1749  areaMeshSubDir,
1750  runTime,
1752  )
1753  );
1754 
1755  // Check faceProcAddressing
1756  labelIOList faceProcAddressing
1757  (
1758  IOobject
1759  (
1760  "faceProcAddressing",
1761  areaMeshInstance,
1762  areaMeshSubDir,
1763  runTime,
1765  )
1766  );
1767 
1768  haveAreaAddressing =
1769  (
1770  faceProcAddressing.headerOk()
1771  && faceProcAddressing.size() == faceLabels.size()
1772  );
1773  }
1774  else if (areaMeshDetected)
1775  {
1776  // Have no mesh. Don't need addressing
1777  haveAreaAddressing = true;
1778  }
1779 
1780 
1781  // Additionally check for master faces being readable. Could
1782  // do even more checks, e.g. global number of cells same
1783  // as cellProcAddressing
1784 
1785  bool volMeshHaveUndecomposed = false;
1786  bool areaMeshHaveUndecomposed = false;
1787 
1788  if (Pstream::master())
1789  {
1790  InfoOrPout
1791  << "Checking " << baseRunTime.caseName()
1792  << " for undecomposed volume and area meshes..."
1793  << endl;
1794 
1795  const bool oldParRun = Pstream::parRun(false);
1796  const label oldNumProcs(fileHandler().nProcs());
1797 
1798  // Volume
1799  {
1800  faceCompactIOList facesIO
1801  (
1802  IOobject
1803  (
1804  "faces",
1805  volMeshInstance,
1806  volMeshSubDir,
1807  baseRunTime,
1809  ),
1810  label(0)
1811  );
1812  volMeshHaveUndecomposed = facesIO.headerOk();
1813  }
1814 
1815  // Area
1816  if (doFiniteArea)
1817  {
1818  labelIOList labelsIO
1819  (
1820  IOobject
1821  (
1822  "faceLabels",
1823  areaMeshInstance,
1824  areaMeshSubDir,
1825  baseRunTime,
1827  )
1828  );
1829  areaMeshHaveUndecomposed = labelsIO.headerOk();
1830  }
1831 
1832  const_cast<fileOperation&>
1833  (
1834  fileHandler()
1835  ).nProcs(oldNumProcs);
1836  Pstream::parRun(oldParRun); // Restore parallel state
1837  }
1838 
1840  (
1842  volMeshHaveUndecomposed,
1843  areaMeshHaveUndecomposed
1844  );
1845 
1846  // Report
1847  {
1848  InfoOrPout
1849  << " volume mesh ["
1850  << volMeshHaveUndecomposed << "] : "
1851  << volMeshInstance << nl
1852  << " area mesh ["
1853  << areaMeshHaveUndecomposed << "] : "
1854  << areaMeshInstance << nl
1855  << endl;
1856  }
1857 
1858 
1859  if
1860  (
1861  !volMeshHaveUndecomposed
1862  || !returnReduceAnd(haveVolAddressing)
1863  )
1864  {
1865  InfoOrPout
1866  << "No undecomposed mesh. Creating from: "
1867  << volMeshInstance << endl;
1868 
1869  if (areaMeshHaveUndecomposed)
1870  {
1871  areaMeshHaveUndecomposed = false;
1872  InfoOrPout
1873  << "Also ignore any undecomposed area mesh"
1874  << endl;
1875  }
1876 
1878  (
1879  IOobject
1880  (
1881  regionName,
1882  volMeshInstance,
1883  runTime,
1885  ),
1886  volMeshReadHandler
1887  );
1888  fvMeshTools::setBasicGeometry(volMeshPtr());
1889  fvMesh& mesh = volMeshPtr();
1890 
1891 
1892  InfoOrPout<< nl << "Reconstructing mesh" << nl << endl;
1893 
1894  // Reconstruct (1 processor)
1895  const label nDestProcs(1);
1896  const labelList finalDecomp(mesh.nCells(), Zero);
1897 
1898  redistributeAndWrite
1899  (
1900  volMeshReadHandler,
1901  masterOnlyHandler, //writeHandler,
1902  baseRunTime,
1903  proc0CaseName,
1904 
1905  // Controls
1906  false, // do not read fields
1907  false, // do not read undecomposed case on proc0
1908  true, // write redistributed files to proc0
1909  overwrite,
1910 
1911  // Decomposition information
1912  nDestProcs,
1913  finalDecomp,
1914 
1915  // For finite-volume
1916  volMeshOnProc,
1917  volMeshInstance,
1918  mesh
1919  );
1920  }
1921 
1922 
1923  // Similarly for finiteArea
1924  // - may or may not have undecomposed mesh
1925  // - may or may not have decomposed meshes
1926 
1927  if
1928  (
1929  areaMeshOnProc.found(true) // ie, areaMeshDetected
1930  &&
1931  (
1932  !areaMeshHaveUndecomposed
1933  || !returnReduceAnd(haveAreaAddressing)
1934  )
1935  )
1936  {
1937  InfoOrPout
1938  << "Loading area mesh from "
1939  << areaMeshInstance << endl;
1940 
1941  InfoOrPout<< " getting volume mesh support" << endl;
1942 
1944  (
1945  IOobject
1946  (
1947  regionName,
1948  baseRunTime.timeName(),
1949  baseRunTime,
1951  ),
1952  true // read on master only
1953  );
1954  fvMeshTools::setBasicGeometry(baseMeshPtr());
1955 
1957  (
1958  IOobject
1959  (
1960  regionName,
1961  baseMeshPtr().facesInstance(),
1962  runTime,
1964  ),
1965  volMeshReadHandler
1966  );
1967  fvMesh& mesh = volMeshPtr();
1968 
1969  // Read volume proc addressing back to base mesh
1971  (
1973  );
1974 
1975 
1976  //autoPtr<faMesh> areaMeshPtr =
1977  // faMeshTools::loadOrCreateMesh
1978  //(
1979  // IOobject
1980  // (
1981  // regionName,
1982  // areaMeshInstance,
1983  // runTime,
1984  // IOobjectOption::MUST_READ
1985  // ),
1986  // mesh, // <- The referenced polyMesh (from above)
1987  // decompose
1988  //);
1990  (
1991  IOobject
1992  (
1993  regionName,
1994  areaMeshInstance,
1995  runTime,
1997  ),
1998  mesh, // <- The referenced polyMesh (from above)
1999  areaMeshReadHandler
2000  );
2001  faMesh& areaMesh = areaMeshPtr();
2002 
2005 
2006  autoPtr<faMesh> areaBaseMeshPtr;
2007 
2008  // Reconstruct using polyMesh distribute map
2009  mapDistributePolyMesh faDistMap
2010  (
2012  (
2013  areaMesh,
2014  distMap(), // The polyMesh distMap
2015  baseMeshPtr(), // Target polyMesh
2016  areaBaseMeshPtr
2017  )
2018  );
2019 
2020  faMeshTools::forceDemandDriven(areaBaseMeshPtr());
2021  faMeshTools::unregisterMesh(areaBaseMeshPtr());
2022 
2023 
2024  if (Pstream::master())
2025  {
2026  InfoOrPout
2027  << "Setting caseName to " << baseRunTime.caseName()
2028  << " to write reconstructed area mesh." << endl;
2029  runTime.caseName() = baseRunTime.caseName();
2030  const bool oldProcCase(runTime.processorCase(false));
2031  const bool oldParRun(Pstream::parRun(false));
2032  const label oldNumProcs(fileHandler().nProcs());
2033 
2034  areaBaseMeshPtr().write();
2035 
2036  // Now we've written all. Reset caseName on master
2037  InfoOrPout<< "Restoring caseName" << endl;
2038  const_cast<fileOperation&>
2039  (
2040  fileHandler()
2041  ).nProcs(oldNumProcs);
2042  Pstream::parRun(oldParRun);
2043  runTime.caseName() = proc0CaseName;
2044  runTime.processorCase(oldProcCase);
2045  }
2046 
2047  // Update for the reconstructed procAddressing
2049  (
2050  areaBaseMeshPtr(), // Reconstruct location
2051  faDistMap,
2052  false, // decompose=false
2053  writeHandler, //writeHandler,
2054  areaMeshPtr.get() // procMesh
2055  );
2056  }
2057  }
2058 
2059  // Make sure all is finished writing until re-reading in pass2
2060  // below
2061  fileHandler().flush();
2062 
2063 
2064  // Pass2 : read mesh and addressing and reconstruct fields
2065  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2066 
2067  InfoOrPout
2068  << nl
2069  << "Reconstructing fields" << nl << endl;
2070 
2071  runTime.setTime(timeDirs[0], 0);
2072  baseRunTime.setTime(timeDirs[0], 0);
2073  InfoOrPout<< "Time = " << runTime.timeName() << endl << endl;
2074 
2075 
2076  // Read undecomposed mesh on master and 'empty' mesh
2077  // (zero faces, point, cells but valid patches and zones) on slaves.
2078  // This is a bit of tricky code and hidden inside fvMeshTools for
2079  // now.
2080  InfoOrPout<< "Reading undecomposed mesh (on master)" << endl;
2081  //autoPtr<fvMesh> baseMeshPtr = fvMeshTools::newMesh
2082  //(
2083  // IOobject
2084  // (
2085  // regionName,
2086  // baseRunTime.timeName(),
2087  // baseRunTime,
2088  // IOobjectOption::MUST_READ
2089  // ),
2090  // true // read on master only
2091  //);
2092  fileName facesInstance;
2093  fileName pointsInstance;
2095  (
2096  IOobject
2097  (
2098  regionName,
2099  baseRunTime.timeName(),
2100  baseRunTime,
2102  ),
2103  facesInstance,
2104  pointsInstance
2105  );
2106 
2108  (
2109  IOobject
2110  (
2111  regionName,
2112  facesInstance, //baseRunTime.timeName(),
2113  baseRunTime,
2115  ),
2116  masterOnlyHandler // read on master only
2117  );
2118 
2119  if (::debug)
2120  {
2121  Pout<< "Undecomposed mesh :"
2122  << " instance:" << baseMeshPtr().facesInstance()
2123  << " nCells:" << baseMeshPtr().nCells()
2124  << " nFaces:" << baseMeshPtr().nFaces()
2125  << " nPoints:" << baseMeshPtr().nPoints()
2126  << endl;
2127  }
2128 
2129  InfoOrPout<< "Reading local, decomposed mesh" << endl;
2131  (
2132  IOobject
2133  (
2134  regionName,
2135  baseMeshPtr().facesInstance(),
2136  runTime,
2138  ),
2139  volMeshReadHandler // read on fvMesh processors
2140  );
2141  fvMesh& mesh = volMeshPtr();
2142 
2143 
2144  // Similarly for finiteArea
2145  autoPtr<faMesh> areaBaseMeshPtr;
2146  autoPtr<faMesh> areaMeshPtr;
2147  autoPtr<faMeshDistributor> faDistributor;
2148  mapDistributePolyMesh areaDistMap;
2149 
2150  if (areaMeshDetected)
2151  {
2152  //areaBaseMeshPtr = faMeshTools::newMesh
2153  //(
2154  // IOobject
2155  // (
2156  // regionName,
2157  // baseRunTime.timeName(),
2158  // baseRunTime,
2159  // IOobjectOption::MUST_READ
2160  // ),
2161  // baseMeshPtr(),
2162  // true // read on master only
2163  //);
2164  areaBaseMeshPtr = faMeshTools::loadOrCreateMesh
2165  (
2166  IOobject
2167  (
2168  regionName,
2169  baseMeshPtr().facesInstance(),
2170  baseRunTime,
2172  ),
2173  baseMeshPtr(),
2174  masterOnlyHandler
2175  );
2176 
2177  //areaMeshPtr = faMeshTools::loadOrCreateMesh
2178  //(
2179  // IOobject
2180  // (
2181  // regionName,
2182  // areaBaseMeshPtr().facesInstance(),
2183  // runTime,
2184  // IOobjectOption::MUST_READ
2185  // ),
2186  // mesh,
2187  // decompose
2188  //);
2189  areaMeshPtr = faMeshTools::loadOrCreateMesh
2190  (
2191  IOobject
2192  (
2193  regionName,
2194  areaBaseMeshPtr().facesInstance(),
2195  runTime,
2197  ),
2198  mesh,
2199  areaMeshReadHandler
2200  );
2201 
2202  areaDistMap =
2204  (
2205  areaMeshPtr(),
2206  areaBaseMeshPtr
2207  );
2208 
2209  faMeshTools::forceDemandDriven(areaMeshPtr());
2210 
2211  // Create an appropriate field distributor
2212  faDistributor.reset
2213  (
2214  new faMeshDistributor
2215  (
2216  areaMeshPtr(), // source
2217  areaBaseMeshPtr(), // target
2218  areaDistMap,
2219  masterOnlyHandler // only write on master
2220  )
2221  );
2222  // Report some messages. Tbd.
2224  }
2225 
2226 
2227  // Read addressing back to base mesh
2229  distMap = fvMeshTools::readProcAddressing(mesh, baseMeshPtr);
2230 
2231  // Construct field mapper
2232  auto fvDistributorPtr =
2234  (
2235  mesh, // source
2236  baseMeshPtr(), // target
2237  distMap(),
2238  masterOnlyHandler // Write on master only
2239  );
2240 
2241  // Construct point field mapper
2242  const auto& basePointMesh = pointMesh::New(baseMeshPtr());
2243  const auto& procPointMesh = pointMesh::New(mesh);
2244 
2245  auto pointFieldDistributorPtr =
2247  (
2248  procPointMesh, // source
2249  basePointMesh, // target
2250  distMap(),
2251  false, // delay
2252  //UPstream::master() // Write reconstructed on master
2253  masterOnlyHandler // Write on master only
2254  );
2255 
2256 
2257  // Since we start from Times[0] and not runTime.timeName() we
2258  // might overlook point motion in the first timestep
2259  // (since mesh.readUpdate() below will not be triggered). Instead
2260  // detect points by hand
2262  {
2263  InfoOrPout
2264  << " Detected initial mesh motion;"
2265  << " reconstructing points" << nl
2266  << endl;
2267  fvDistributorPtr().reconstructPoints();
2268  }
2269 
2270 
2271  // Loop over all times
2272  forAll(timeDirs, timeI)
2273  {
2274  if (newTimes && masterTimeDirSet.found(timeDirs[timeI].name()))
2275  {
2276  InfoOrPout
2277  << "Skipping time " << timeDirs[timeI].name()
2278  << nl << endl;
2279  continue;
2280  }
2281 
2282  // Set time for global database
2283  runTime.setTime(timeDirs[timeI], timeI);
2284  baseRunTime.setTime(timeDirs[timeI], timeI);
2285 
2286  InfoOrPout<< "Time = " << runTime.timeName() << endl << endl;
2287 
2288 
2289  // Check if any new meshes need to be read.
2291 
2292  if (procStat == polyMesh::POINTS_MOVED)
2293  {
2294  InfoOrPout
2295  << " Detected mesh motion; reconstructing points"
2296  << nl << endl;
2297  fvDistributorPtr().reconstructPoints();
2298  }
2299  else if
2300  (
2301  procStat == polyMesh::TOPO_CHANGE
2302  || procStat == polyMesh::TOPO_PATCH_CHANGE
2303  )
2304  {
2305  InfoOrPout
2306  << " Detected topology change;"
2307  << " reconstructing addressing" << nl << endl;
2308 
2309  if (baseMeshPtr)
2310  {
2311  // Cannot do a baseMesh::readUpdate() since not all
2312  // processors will have mesh files. So instead just
2313  // recreate baseMesh
2314  baseMeshPtr.clear();
2315  //baseMeshPtr = fvMeshTools::newMesh
2316  //(
2317  // IOobject
2318  // (
2319  // regionName,
2320  // baseRunTime.timeName(),
2321  // baseRunTime,
2322  // IOobjectOption::MUST_READ
2323  // ),
2324  // true // read on master only
2325  //);
2326  baseMeshPtr = fvMeshTools::loadOrCreateMesh
2327  (
2328  IOobject
2329  (
2330  regionName,
2331  baseRunTime.timeName(),
2332  baseRunTime,
2334  ),
2335  masterOnlyHandler // read on master only
2336  );
2337  if (::debug)
2338  {
2339  Pout<< "Undecomposed mesh :"
2340  << " nCells:" << baseMeshPtr().nCells()
2341  << " nFaces:" << baseMeshPtr().nFaces()
2342  << " nPoints:" << baseMeshPtr().nPoints()
2343  << endl;
2344  }
2345  }
2346 
2347  // Re-read procAddressing
2348  distMap =
2349  fvMeshTools::readProcAddressing(mesh, baseMeshPtr);
2350 
2351  // Reset field mappers
2352 
2353  fvDistributorPtr.reset
2354  (
2356  (
2357  mesh, // source
2358  baseMeshPtr(), // target
2359  distMap(),
2360  masterOnlyHandler // Write on master only
2361  )
2362  );
2363 
2364  // Construct point field mapper
2365  const auto& basePointMesh = pointMesh::New(baseMeshPtr());
2366  const auto& procPointMesh = pointMesh::New(mesh);
2367 
2368  pointFieldDistributorPtr.reset
2369  (
2371  (
2372  procPointMesh, // source
2373  basePointMesh, // target
2374  distMap(),
2375  false, // delay until later
2376  masterOnlyHandler // Write on master only
2377  )
2378  );
2379 
2380  lagrangianDistributorPtr.reset();
2381 
2382  if (areaMeshPtr)
2383  {
2384  InfoOrPout
2385  << " Discarding finite-area addressing"
2386  << " (TODO)" << nl << endl;
2387 
2388  areaBaseMeshPtr.reset();
2389  areaMeshPtr.reset();
2390  faDistributor.reset();
2391  areaDistMap.clear();
2392  }
2393  }
2394 
2395 
2396  // Get list of objects
2397  IOobjectList objects(mesh, runTime.timeName());
2398 
2399  // Mesh fields (vol, surface, volInternal)
2400  fvDistributorPtr()
2401  .distributeAllFields(objects, selectedFields);
2402 
2403  // pointfields
2404  // - distribute and write (verbose)
2405  pointFieldDistributorPtr()
2406  .distributeAllFields(objects, selectedFields);
2407 
2408 
2409  // Clouds (note: might not be present on all processors)
2411  (
2412  lagrangianDistributorPtr,
2413  baseMeshPtr(),
2414  mesh,
2415  distMap(),
2416  selectedLagrangianFields
2417  //masterOnlyHandler
2418  );
2419 
2420  if (faDistributor)
2421  {
2422  faDistributor()
2423  .distributeAllFields(objects, selectedFields);
2424  }
2425 
2426 
2427  // If there are any "uniform" directories copy them from
2428  // the master processor
2429 
2430  copyUniform
2431  (
2432  volMeshReadHandler, //masterOnlyHandler, // readHandler
2433  masterOnlyHandler, // writeHandler
2434 
2435  true, // reconstruct
2436  false, // decompose
2437 
2438  mesh.time().timeName(),
2439  word::null, // optional caseName for reading
2440  mesh,
2441  baseMeshPtr()
2442  );
2443  // Non-region specific. Note: should do outside region loop
2444  // but would then have to replicate the whole time loop ...
2445  copyUniform
2446  (
2447  volMeshReadHandler, //masterOnlyHandler, // readHandler,
2448  masterOnlyHandler, // writeHandler
2449 
2450  true, // reconstruct
2451  false, // decompose
2452 
2453  mesh.time().timeName(),
2454  word::null, // optional caseName for reading
2455  mesh.time(), // runTime
2456  baseMeshPtr().time() // baseRunTime
2457  );
2458  }
2459  }
2460  }
2461  else
2462  {
2463  // decompose or redistribution mode.
2464  // decompose : master : read from parent dir
2465  // slave : dummy mesh
2466  // redistribute : all read mesh or dummy mesh
2467 
2468  Time& readRunTime =
2469  (
2470  (decompose)
2471  ? baseRunTime
2472  : runTime
2473  );
2474 
2475  // Time coming from processor0 (or undecomposed if no processor0)
2476  scalar masterTime = timeSelector::selectIfPresent
2477  (
2478  readRunTime,
2479  args
2480  )[0].value();
2481  Pstream::broadcast(masterTime);
2482  InfoOrPout
2483  << "Setting time to that of master or undecomposed case : "
2484  << masterTime << endl;
2485  runTime.setTime(masterTime, 0);
2486  baseRunTime.setTime(masterTime, 0);
2487 
2488 
2489 
2490 
2491  // Save old time name (since might be incremented)
2492  const word oldTimeName(runTime.timeName());
2493 
2494  forAll(regionNames, regioni)
2495  {
2496  const word& regionName = regionNames[regioni];
2497 
2498  const fileName volMeshSubDir
2499  (
2501  );
2502  const fileName areaMeshSubDir
2503  (
2505  );
2506 
2507  InfoOrPout
2508  << nl << nl
2509  << (decompose ? "Decomposing" : "Redistributing")
2510  << " mesh:" << polyMesh::regionName(regionName) << nl << endl;
2511 
2512 
2513  // Get time instance directory
2514  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
2515  // At this point we should be able to read at least a mesh on
2516  // processor0. Note the changing of the processor0 casename to
2517  // enforce it to read/write from the undecomposed case
2518 
2519  fileName volMeshMasterInstance;
2520  fileName areaMeshMasterInstance;
2521 
2522  // Assume to be true
2523  bool volMeshHaveUndecomposed = true;
2524  bool areaMeshHaveUndecomposed = doFiniteArea;
2525 
2526  if (Pstream::master())
2527  {
2528  if (decompose)
2529  {
2530  InfoOrPout
2531  << "Checking undecomposed mesh in case: "
2532  << baseRunTime.caseName() << endl;
2533  runTime.caseName() = baseRunTime.caseName();
2534  runTime.processorCase(false);
2535  }
2536 
2537  const bool oldParRun = Pstream::parRun(false);
2538  const label oldNumProcs(fileHandler().nProcs());
2539  volMeshMasterInstance = readRunTime.findInstance
2540  (
2541  volMeshSubDir,
2542  "faces",
2544  );
2545 
2546  if (doFiniteArea)
2547  {
2548  areaMeshMasterInstance = readRunTime.findInstance
2549  (
2550  areaMeshSubDir,
2551  "faceLabels",
2553  );
2554 
2555  // Note: findInstance returns "constant" even if not found,
2556  // so recheck now for a false positive.
2557 
2558  if ("constant" == areaMeshMasterInstance)
2559  {
2560  const boolList areaMeshOnProc
2561  (
2562  haveMeshFile
2563  (
2564  readRunTime,
2565  areaMeshMasterInstance/areaMeshSubDir,
2566  "faceLabels",
2567  false // verbose=false
2568  )
2569  );
2570 
2571  if (areaMeshOnProc.empty() || !areaMeshOnProc[0])
2572  {
2573  areaMeshHaveUndecomposed = false;
2574  }
2575  }
2576  }
2577 
2578  const_cast<fileOperation&>(fileHandler()).nProcs(oldNumProcs);
2579  Pstream::parRun(oldParRun); // Restore parallel state
2580 
2581  if (decompose)
2582  {
2583  InfoOrPout
2584  << " volume mesh ["
2585  << volMeshHaveUndecomposed << "] : "
2586  << volMeshMasterInstance << nl
2587  << " area mesh ["
2588  << areaMeshHaveUndecomposed << "] : "
2589  << areaMeshMasterInstance << nl
2590  << nl << nl;
2591 
2592  // Restoring caseName
2593  InfoOrPout<< "Restoring caseName" << endl;
2594  runTime.caseName() = proc0CaseName;
2595  runTime.processorCase(oldProcCase);
2596  }
2597  }
2598 
2600  (
2602  volMeshHaveUndecomposed,
2603  areaMeshHaveUndecomposed,
2604  volMeshMasterInstance,
2605  areaMeshMasterInstance
2606  );
2607 
2608  // Check processors have meshes
2609  // - check for 'faces' file (polyMesh)
2610  // - check for 'faceLabels' file (faMesh)
2611  boolList volMeshOnProc;
2612  boolList areaMeshOnProc;
2613 
2614  if (decompose)
2615  {
2616  // Already determined above that master can read 'faces' file.
2617  // This avoids doing all the casename setting/restoring again.
2618  volMeshOnProc.setSize(UPstream::nProcs(), false);
2619  volMeshOnProc[UPstream::masterNo()] = volMeshHaveUndecomposed;
2620  }
2621  else
2622  {
2623  // All check if can read 'faces' file
2624  volMeshOnProc = haveMeshFile
2625  (
2626  runTime,
2627  volMeshMasterInstance/volMeshSubDir,
2628  "faces"
2629  );
2630  }
2631 
2632  // Create handler for reading
2633  newHandler(volMeshOnProc, volMeshReadHandler);
2634 
2635 
2636  // Now we've determined which processors are reading switch back
2637  // to exact matching of 'processorsXXX' directory names.
2638  // - this determines the processorsXXX fileNames
2639  // - the XXX comes from the number of read processors
2640  // - also adapt the masterOnlyReader (used in copyUniform)
2641 
2643  (
2644  findIndices(volMeshOnProc, true).size()
2645  );
2646 
2647 
2648  if (doFiniteArea)
2649  {
2650  if (decompose)
2651  {
2652  // Already determined above that master can read
2653  // 'faceLabels' file.
2654  areaMeshOnProc.setSize(UPstream::nProcs(), false);
2655  areaMeshOnProc[UPstream::masterNo()] =
2656  areaMeshHaveUndecomposed;
2657  }
2658  else
2659  {
2660  areaMeshOnProc = haveMeshFile
2661  (
2662  runTime,
2663  areaMeshMasterInstance/areaMeshSubDir,
2664  "faceLabels"
2665  );
2666  }
2667 
2668  // Create handler for reading
2669  if (areaMeshOnProc == volMeshOnProc)
2670  {
2671  if (volMeshReadHandler)
2672  {
2673  // Use same reader for faMesh as for fvMesh
2674  areaMeshReadHandler.ref(volMeshReadHandler.ref());
2675  }
2676  }
2677  else
2678  {
2679  newHandler(areaMeshOnProc, areaMeshReadHandler);
2680  }
2681  }
2682 
2683 
2684  // Prior to loadOrCreateMesh, note which meshes already exist
2685  // for the current file handler.
2686  // - where mesh would be written if it didn't exist already.
2687  fileNameList volMeshDir(Pstream::nProcs());
2688  {
2689  volMeshDir[Pstream::myProcNo()] =
2690  (
2691  fileHandler().objectPath
2692  (
2693  IOobject
2694  (
2695  "faces",
2696  volMeshMasterInstance/volMeshSubDir,
2697  runTime
2698  ),
2699  word::null
2700  ).path()
2701  );
2702 
2703  Pstream::allGatherList(volMeshDir);
2704 
2705  if (optVerbose && Pstream::master())
2706  {
2707  Info<< "Per processor faces dirs:" << nl
2708  << '(' << nl;
2709 
2710  for (const int proci : Pstream::allProcs())
2711  {
2712  Info<< " "
2713  << runTime.relativePath(volMeshDir[proci]);
2714 
2715  if (!volMeshOnProc[proci])
2716  {
2717  Info<< " [missing]";
2718  }
2719  Info<< nl;
2720  }
2721  Info<< ')' << nl << endl;
2722  }
2723  }
2724 
2725  fileNameList areaMeshDir(Pstream::nProcs());
2726  if (doFiniteArea)
2727  {
2728  areaMeshDir[Pstream::myProcNo()] =
2729  (
2730  fileHandler().objectPath
2731  (
2732  IOobject
2733  (
2734  "faceLabels",
2735  areaMeshMasterInstance/areaMeshSubDir,
2736  runTime
2737  ),
2738  word::null
2739  ).path()
2740  );
2741 
2742  Pstream::allGatherList(areaMeshDir);
2743 
2744  if (optVerbose && Pstream::master())
2745  {
2746  Info<< "Per processor faceLabels dirs:" << nl
2747  << '(' << nl;
2748 
2749  for (const int proci : Pstream::allProcs())
2750  {
2751  Info<< " "
2752  << runTime.relativePath(areaMeshDir[proci]);
2753 
2754  if (!areaMeshOnProc[proci])
2755  {
2756  Info<< " [missing]";
2757  }
2758  Info<< nl;
2759  }
2760  Info<< ')' << nl << endl;
2761  }
2762  }
2763 
2764 
2765  // Load mesh (or create dummy one)
2766  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2767 
2768  if (decompose)
2769  {
2770  InfoOrPout
2771  << "Setting caseName to " << baseRunTime.caseName()
2772  << " to read undecomposed mesh" << endl;
2773  runTime.caseName() = baseRunTime.caseName();
2774  runTime.processorCase(false);
2775  }
2776 
2778  (
2779  IOobject
2780  (
2781  regionName,
2782  volMeshMasterInstance,
2783  runTime,
2785  ),
2786  volMeshReadHandler
2787  );
2788  fvMesh& mesh = volMeshPtr();
2789 
2790 
2791  // Area mesh
2792 
2793  autoPtr<faMesh> areaMeshPtr;
2794 
2795  // Decomposing: must have an undecomposed mesh
2796  // Redistributing: have any proc mesh
2797  if
2798  (
2799  doFiniteArea
2800  &&
2801  (
2802  decompose
2803  ? areaMeshHaveUndecomposed
2804  : areaMeshOnProc.found(true)
2805  )
2806  )
2807  {
2808  areaMeshPtr = faMeshTools::loadOrCreateMesh
2809  (
2810  IOobject
2811  (
2812  regionName,
2813  areaMeshMasterInstance,
2814  runTime,
2816  ),
2817  mesh, // <- The referenced polyMesh (from above)
2818  areaMeshReadHandler
2819  );
2820 
2821  faMeshTools::forceDemandDriven(*areaMeshPtr);
2822  faMeshTools::unregisterMesh(*areaMeshPtr);
2823  }
2824 
2825 
2826  if (decompose)
2827  {
2828  InfoOrPout<< "Restoring caseName" << endl;
2829  runTime.caseName() = proc0CaseName;
2830  runTime.processorCase(oldProcCase);
2831  }
2832 
2833  const label nOldCells = mesh.nCells();
2834 
2835 
2836  // Determine decomposition
2837  // ~~~~~~~~~~~~~~~~~~~~~~~
2838 
2839  label nDestProcs;
2840  labelList finalDecomp;
2841  determineDecomposition
2842  (
2843  volMeshReadHandler, // how to read decomposeParDict
2844  baseRunTime,
2845  decompDictFile,
2846  decompose,
2847  proc0CaseName,
2848  mesh,
2849  writeCellDist,
2850 
2851  nDestProcs,
2852  finalDecomp
2853  );
2854 
2855  if (dryrun)
2856  {
2857  continue;
2858  }
2859 
2860  if (!writeHandler && nDestProcs < fileHandler().nProcs())
2861  {
2862  boolList isWriteProc(UPstream::nProcs(), false);
2863  isWriteProc.slice(0, nDestProcs) = true;
2864  InfoOrPout
2865  << " dest procs ["
2866  << isWriteProc << "]" << nl
2867  << endl;
2868  newHandler(isWriteProc, writeHandler);
2869  }
2870 
2871  // Area fields first. Read and deregister
2873  if (areaMeshPtr)
2874  {
2875  areaFields.read
2876  (
2877  baseRunTime,
2878  proc0CaseName,
2879  decompose,
2880 
2881  areaMeshOnProc,
2882  areaMeshReadHandler,
2883  areaMeshMasterInstance,
2884  (*areaMeshPtr)
2885  );
2886  }
2887 
2888 
2889  // Detect lagrangian fields
2890  if (decompose)
2891  {
2892  InfoOrPout
2893  << "Setting caseName to " << baseRunTime.caseName()
2894  << " to read lagrangian" << endl;
2895  if (UPstream::master())
2896  {
2897  // Change case name but only on the master - this will
2898  // hopefully cause the slaves to not read.
2899  runTime.caseName() = baseRunTime.caseName();
2900  }
2901  else
2902  {
2903  // Explicitly make sure that casename is not recognised as
2904  // a processor case since that has special handling for
2905  // caching processor directories etc.
2906  runTime.caseName() = "#invalid-name#";
2907  }
2908  runTime.processorCase(false);
2909  }
2910 
2911  // Read lagrangian fields and store on cloud (objectRegistry)
2913  (
2915  (
2916  mesh,
2917  selectedLagrangianFields
2918  )
2919  );
2920 
2921  if (decompose)
2922  {
2923  InfoOrPout<< "Restoring caseName" << endl;
2924  runTime.caseName() = proc0CaseName;
2925  runTime.processorCase(oldProcCase);
2926  }
2927 
2928 
2929  // Load fields, do all distribution (mesh and fields)
2930  // - but not lagrangian fields; these are done later
2931  autoPtr<mapDistributePolyMesh> distMap = redistributeAndWrite
2932  (
2933  volMeshReadHandler, // readHandler
2934  writeHandler, // writeHandler,
2935  baseRunTime,
2936  proc0CaseName,
2937 
2938  // Controls
2939  true, // read fields
2940  decompose, // decompose, i.e. read from undecomposed case
2941  false, // no reconstruction
2942  overwrite,
2943 
2944  // Decomposition information
2945  nDestProcs,
2946  finalDecomp,
2947 
2948  // For finite volume
2949  volMeshOnProc,
2950  volMeshMasterInstance,
2951  mesh
2952  );
2953 
2954  // Redistribute any clouds
2956  (
2957  lagrangianDistributorPtr,
2958  mesh,
2959  nOldCells,
2960  distMap(),
2961  clouds
2962  );
2963 
2964 
2965  // Redistribute area fields
2966 
2967  mapDistributePolyMesh faDistMap;
2968  autoPtr<faMesh> areaProcMeshPtr;
2969 
2970  if (areaMeshPtr)
2971  {
2972  faDistMap = faMeshDistributor::distribute
2973  (
2974  areaMeshPtr(),
2975  distMap(),
2976  areaProcMeshPtr
2977  );
2978 
2979  // Force recreation of everything that might vaguely
2980  // be used by patches:
2981 
2982  faMeshTools::forceDemandDriven(areaProcMeshPtr());
2983 
2984 
2985  if (reconstruct)
2986  {
2987  if (Pstream::master())
2988  {
2989  InfoOrPout
2990  << "Setting caseName to " << baseRunTime.caseName()
2991  << " to write reconstructed mesh (and fields)."
2992  << endl;
2993  runTime.caseName() = baseRunTime.caseName();
2994  const bool oldProcCase(runTime.processorCase(false));
2995  const bool oldParRun = UPstream::parRun(false);
2996  const label oldNumProcs(fileHandler().nProcs());
2997 
2998  areaProcMeshPtr->write();
2999 
3000  // Now we've written all. Reset caseName on master
3001  InfoOrPout<< "Restoring caseName" << endl;
3002  const_cast<fileOperation&>
3003  (
3004  fileHandler()
3005  ).nProcs(oldNumProcs);
3006  UPstream::parRun(oldParRun);
3007  runTime.caseName() = proc0CaseName;
3008  runTime.processorCase(oldProcCase);
3009  }
3010  }
3011  else
3012  {
3013  auto oldHandler = fileOperation::fileHandler(writeHandler);
3014 
3016  (
3017  IOobject
3018  (
3019  "procAddressing",
3020  areaProcMeshPtr->facesInstance(),
3022  areaProcMeshPtr->thisDb(),
3026  ),
3027  faDistMap
3028  ).write();
3029 
3030  areaProcMeshPtr->write();
3031 
3032  writeHandler = fileOperation::fileHandler(oldHandler);
3033 
3034  if (decompose)
3035  {
3037  (
3038  areaProcMeshPtr(),
3039  faDistMap,
3040  decompose,
3041  writeHandler
3042  );
3043  }
3044  }
3045 
3046  InfoOrPout
3047  << "Written redistributed mesh to "
3048  << areaProcMeshPtr->facesInstance() << nl << endl;
3049 
3050  faMeshDistributor distributor
3051  (
3052  areaMeshPtr(), // source
3053  areaProcMeshPtr(), // target
3054  faDistMap,
3055  writeHandler
3056  );
3057 
3058  areaFields.redistributeAndWrite(distributor, true);
3059  }
3060 
3061 
3062  // Get reference to standard write handler
3063  refPtr<fileOperation> defaultHandler;
3064  if (writeHandler)
3065  {
3066  defaultHandler.ref(writeHandler.ref());
3067  }
3068  else
3069  {
3070  defaultHandler.ref(const_cast<fileOperation&>(fileHandler()));
3071  }
3072 
3073 
3074  copyUniform
3075  (
3076  volMeshReadHandler, // read handler
3077  defaultHandler, //TBD: should be all IOranks
3078 
3079  reconstruct, // reconstruct
3080  decompose, // decompose
3081 
3082  oldTimeName,
3083  (decompose ? baseRunTime.caseName() : proc0CaseName),
3084  mesh, // read location is mesh (but oldTimeName)
3085  mesh // write location is mesh
3086  );
3087  }
3088 
3089 
3090  // Get reference to standard write handler
3091  refPtr<fileOperation> defaultHandler;
3092  if (writeHandler)
3093  {
3094  defaultHandler.ref(writeHandler.ref());
3095  }
3096  else
3097  {
3098  defaultHandler.ref(const_cast<fileOperation&>(fileHandler()));
3099  }
3100 
3101  copyUniform
3102  (
3103  volMeshReadHandler, // read handler
3104  defaultHandler, //TBD: should be all IOranks
3105 
3106  reconstruct, // reconstruct (=false)
3107  decompose, // decompose
3108 
3109  oldTimeName, // provided read time
3110  (decompose ? baseRunTime.caseName() : proc0CaseName),
3111  readRunTime,
3112  runTime // writing location
3113  );
3114  }
3115 
3116 
3117  InfoOrPout<< "End\n" << endl;
3118 
3119  return 0;
3120 }
3121 
3122 
3123 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
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:87
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
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")
A class for handling file names.
Definition: fileName.H:72
readOption readOpt() const noexcept
Get the read option.
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
static int nProcsFilter() noexcept
Return collated &#39;processorsDDD&#39; filtering.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
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:725
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:598
bool found(const Key &key) const
Same as contains()
Definition: HashTable.H:1354
static void unset(bool verbose=false)
Deactivate SIGFPE handler and NaN memory initialisation.
Definition: sigFpe.C:208
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:847
static void writeProcAddressing(const fvMesh &procMesh, const mapDistributePolyMesh &map, const bool decompose, const fileName &writeHandlerInstance, refPtr< fileOperation > &writeHandler)
Write addressing if decomposing (1 to many) or reconstructing (many to 1)
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:1176
const fileName & caseName() const noexcept
Return case name.
Definition: TimePathsI.H:78
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
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:666
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:410
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool found(const T &val, label pos=0) const
Same as contains()
Definition: UList.H:879
wordList names() const
The unsorted names of the IOobjects.
Definition: IOobjectList.C:218
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
void clear() noexcept
Same as reset(nullptr)
Definition: autoPtr.H:255
bool erase(iterator &iter)
Erase entry specified by given iterator and delete the allocated pointer.
Definition: HashPtrTable.C:109
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
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 the rootPath.
Definition: TimePathsI.H:66
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...
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
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:1074
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:360
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:409
A IOmapDistributePolyMesh wrapper for using referenced external data.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
labelList faceLabels(nFaceLabels)
Distributor/redistributor for point fields, uses a two (or three) stage construction.
static word processorsBaseDir
Return the processors directory name (usually "processors")
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:815
bool processorCase() const noexcept
True if this is a processor case.
Definition: TimePathsI.H:52
An encapsulation of filesystem-related operations.
static mapDistributePolyMesh readProcAddressing(const faMesh &mesh, const autoPtr< faMesh > &baseMeshPtr)
Read decompose/reconstruct addressing.
UPtrList< Type > sorted()
Return sorted list of objects with a class satisfying isA<Type> or isType<Type> (with Strict) ...
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
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
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
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:421
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. ...
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:61
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:1065
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:853
void setSize(const label n)
Alias for resize()
Definition: List.H:316
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 expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:704
void masterMeshInstance(const IOobject &io, fileName &facesInstance, fileName &pointsInstance)
Determine master faces instance.
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 functionObjectList & functionObjects() const noexcept
Return the list of function objects.
Definition: Time.H:722
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
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:405
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:118
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition: UPstream.H:1059
const word & system() const noexcept
Return system name.
Definition: TimePathsI.H:118
static int verbose_
Output verbosity when writing.
Reading is optional [identical to LAZY_READ].
Simple container to manage read/write, redistribute finiteArea fields.
const fileName & globalCaseName() const noexcept
Return global case name.
Definition: TimePathsI.H:72
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
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
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
Finite area area (element) fields.
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
label localSize(const label proci) const
Size of proci data.
Definition: globalIndexI.H:257
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:935
Abstract base class for domain decomposition.
MeshObject wrapper of decompositionMethod.
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:268
const dictionary & controlDict() const noexcept
Return read access to the controlDict dictionary.
Definition: Time.H:539
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:1177
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
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1102
void reconstructLagrangian(autoPtr< parLagrangianDistributor > &distributorPtr, const fvMesh &baseMesh, const fvMesh &mesh, const mapDistributePolyMesh &distMap, const wordRes &selectedFields)
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:1143
static int cacheLevel() noexcept
Return cache level.
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
T * get() noexcept
Return pointer without nullptr checking.
Definition: refPtr.H:249
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:112
int debug
Static debugging option.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: refPtrI.H:230
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< fileOperation > New(const word &handlerType, bool verbose=false)
Select fileHandler-type. Uses defaultFileHandler if the handlerType is empty.
static label commWorld() noexcept
Communicator for all ranks (respecting any local worlds)
Definition: UPstream.H:429
static autoPtr< faMesh > loadOrCreateMesh(const IOobject &io, const polyMesh &pMesh, const bool decompose, const bool verbose=false)
Definition: faMeshTools.C:580
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional &#39;FOAM Warning&#39; header text...
static void broadcasts(const label comm, Type &arg1, Args &&... args)
Broadcast multiple items to all communicator ranks. Does nothing in non-parallel. ...
static instantList findTimes(const fileName &directory, const word &constantDirName="constant")
Search a given directory for valid time directories.
Definition: TimePaths.C:109
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))
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:701
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: faMesh.C:803
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
static void unregisterMesh(const faMesh &mesh)
Unregister the faMesh from its associated polyMesh to prevent triggering on polyMesh changes etc...
Definition: faMeshTools.C:33
instantList times() const
Search the case for valid time directories.
Definition: TimePaths.C:118
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
#define WarningInFunction
Report a warning using Foam::Warning.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
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:78
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:841
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
Mesh data needed to do the Finite Area discretisation.
Definition: areaFaMesh.H:47
label nPointFields
Required Classes.
Nothing to be read.
Automatically write from objectRegistry::writeObject()
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
Reading is optional [identical to READ_IF_PRESENT].
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:74
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
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
Switch runTimeModifiable() const noexcept
Supports re-reading.
Definition: Time.H:585
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Foam::argList args(argc, argv)
A IOList wrapper for writing external data.
Definition: IOList.H:151
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
static autoPtr< mapDistributePolyMesh > readProcAddressing(const fvMesh &procMesh, const autoPtr< fvMesh > &baseMeshPtr)
Read procAddressing components (reconstructing)
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
PtrList< unmappedPassivePositionParticleCloud > readLagrangian(const fvMesh &mesh, const wordList &cloudNames, const boolUList &haveClouds, const wordRes &selectedFields)
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
virtual labelList decompose(const pointField &points, const scalarField &pointWeights=scalarField::null()) const
Return the wanted processor number for every coordinate, using uniform or specified point weights...
static autoPtr< fileOperation > NewUncollated()
The commonly used uncollatedFileOperation.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127