redistributePar.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  redistributePar
26 
27 Description
28  Redistributes existing decomposed mesh and fields according to the current
29  settings in the decomposeParDict file.
30 
31  Must be run on maximum number of source and destination processors.
32  Balances mesh and writes new mesh to new time directory.
33 
34  Can also work like decomposePar:
35  \verbatim
36  # Create empty processor directories (have to exist for argList)
37  mkdir processor0
38  ..
39  mkdir processorN
40 
41  # Copy undecomposed polyMesh
42  cp -R constant processor0
43 
44  # Distribute
45  mpirun -np ddd redistributePar -parallel
46  \endverbatim
47 \*---------------------------------------------------------------------------*/
48 
49 #include "fvMesh.H"
50 #include "decompositionMethod.H"
51 #include "PstreamReduceOps.H"
52 #include "fvCFD.H"
53 #include "fvMeshDistribute.H"
54 #include "polyDistributionMap.H"
55 #include "IOobjectList.H"
56 #include "globalIndex.H"
57 #include "loadOrCreateMesh.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 void printMeshData(const polyMesh& mesh)
62 {
63  // Collect all data on master
64 
65  globalIndex globalCells(mesh.nCells());
66  labelListList patchNeiProcNo(Pstream::nProcs());
67  labelListList patchSize(Pstream::nProcs());
68  const labelList& pPatches = mesh.globalData().processorPatches();
69  patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
70  patchSize[Pstream::myProcNo()].setSize(pPatches.size());
71  forAll(pPatches, i)
72  {
73  const processorPolyPatch& ppp = refCast<const processorPolyPatch>
74  (
75  mesh.boundaryMesh()[pPatches[i]]
76  );
77  patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
78  patchSize[Pstream::myProcNo()][i] = ppp.size();
79  }
80  Pstream::gatherList(patchNeiProcNo);
81  Pstream::gatherList(patchSize);
82 
83 
84  // Print stats
85 
86  globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());
87 
88  label maxProcCells = 0;
89  label totProcFaces = 0;
90  label maxProcPatches = 0;
91  label totProcPatches = 0;
92  label maxProcFaces = 0;
93 
94  for (label proci = 0; proci < Pstream::nProcs(); proci++)
95  {
96  Info<< endl
97  << "Processor " << proci << nl
98  << " Number of cells = " << globalCells.localSize(proci)
99  << endl;
100 
101  label nProcFaces = 0;
102 
103  const labelList& nei = patchNeiProcNo[proci];
104 
105  forAll(patchNeiProcNo[proci], i)
106  {
107  Info<< " Number of faces shared with processor "
108  << patchNeiProcNo[proci][i] << " = " << patchSize[proci][i]
109  << endl;
110 
111  nProcFaces += patchSize[proci][i];
112  }
113 
114  Info<< " Number of processor patches = " << nei.size() << nl
115  << " Number of processor faces = " << nProcFaces << nl
116  << " Number of boundary faces = "
117  << globalBoundaryFaces.localSize(proci) << endl;
118 
119  maxProcCells = max(maxProcCells, globalCells.localSize(proci));
120  totProcFaces += nProcFaces;
121  totProcPatches += nei.size();
122  maxProcPatches = max(maxProcPatches, nei.size());
123  maxProcFaces = max(maxProcFaces, nProcFaces);
124  }
125 
126  // Stats
127 
128  scalar avgProcCells = scalar(globalCells.size())/Pstream::nProcs();
129  scalar avgProcPatches = scalar(totProcPatches)/Pstream::nProcs();
130  scalar avgProcFaces = scalar(totProcFaces)/Pstream::nProcs();
131 
132  // In case of all faces on one processor. Just to avoid division by 0.
133  if (totProcPatches == 0)
134  {
135  avgProcPatches = 1;
136  }
137  if (totProcFaces == 0)
138  {
139  avgProcFaces = 1;
140  }
141 
142  Info<< nl
143  << "Number of processor faces = " << totProcFaces/2 << nl
144  << "Max number of cells = " << maxProcCells
145  << " (" << 100.0*(maxProcCells-avgProcCells)/avgProcCells
146  << "% above average " << avgProcCells << ")" << nl
147  << "Max number of processor patches = " << maxProcPatches
148  << " (" << 100.0*(maxProcPatches-avgProcPatches)/avgProcPatches
149  << "% above average " << avgProcPatches << ")" << nl
150  << "Max number of faces between processors = " << maxProcFaces
151  << " (" << 100.0*(maxProcFaces-avgProcFaces)/avgProcFaces
152  << "% above average " << avgProcFaces << ")" << nl
153  << endl;
154 }
155 
156 
157 // Debugging: write volScalarField with decomposition for post processing.
158 void writeDecomposition
159 (
160  const word& name,
161  const fvMesh& mesh,
162  const labelList& decomp
163 )
164 {
165  Info<< "Writing wanted cell distribution to volScalarField " << name
166  << " for postprocessing purposes." << nl << endl;
167 
168  volScalarField procCells
169  (
170  IOobject
171  (
172  name,
173  mesh.time().timeName(),
174  mesh,
175  IOobject::NO_READ,
176  IOobject::AUTO_WRITE,
177  false // do not register
178  ),
179  mesh,
180  dimensionedScalar(name, dimless, -1),
181  extrapolatedCalculatedFvPatchScalarField::typeName
182  );
183 
184  forAll(procCells, cI)
185  {
186  procCells[cI] = decomp[cI];
187  }
188  procCells.correctBoundaryConditions();
189 
190  procCells.write();
191 }
192 
193 
194 template<class GeoField>
195 void readFields
196 (
197  const boolList& haveMesh,
198  const typename GeoField::Mesh& mesh,
199  const autoPtr<fvMeshSubset>& subsetterPtr,
200  IOobjectList& allObjects,
201  PtrList<GeoField>& fields
202 )
203 {
204  // Get my objects of type
205  IOobjectList objects(allObjects.lookupClass(GeoField::typeName));
206 
207  // Check that we all have all objects
208  wordList objectNames = objects.toc();
209 
210  // Get master names
211  wordList masterNames(objectNames);
212  Pstream::scatter(masterNames);
213 
214  if (haveMesh[Pstream::myProcNo()] && objectNames != masterNames)
215  {
217  << "differing fields of type " << GeoField::typeName
218  << " on processors." << endl
219  << "Master has:" << masterNames << endl
220  << Pstream::myProcNo() << " has:" << objectNames
221  << abort(FatalError);
222  }
223 
224  fields.setSize(masterNames.size());
225 
226  // Have master send all fields to processors that don't have a mesh
227  if (Pstream::master())
228  {
229  forAll(masterNames, i)
230  {
231  const word& name = masterNames[i];
232  IOobject& io = *objects[name];
233  io.writeOpt() = IOobject::AUTO_WRITE;
234 
235  // Load field
236  fields.set(i, new GeoField(io, mesh));
237 
238  // Create zero sized field and send
239  if (subsetterPtr.valid())
240  {
241  tmp<GeoField> tsubfld = subsetterPtr().interpolate(fields[i]);
242 
243  // Send to all processors that don't have a mesh
244  for (label proci = 1; proci < Pstream::nProcs(); proci++)
245  {
246  if (!haveMesh[proci])
247  {
248  OPstream toProc(Pstream::commsTypes::blocking, proci);
249  toProc<< tsubfld();
250  }
251  }
252  }
253  }
254  }
255  else if (!haveMesh[Pstream::myProcNo()])
256  {
257  // Don't have mesh (nor fields). Receive empty field from master.
258 
259  forAll(masterNames, i)
260  {
261  const word& name = masterNames[i];
262 
263  // Receive field
264  IPstream fromMaster
265  (
266  Pstream::commsTypes::blocking,
267  Pstream::masterNo()
268  );
269  dictionary fieldDict(fromMaster);
270 
271  fields.set
272  (
273  i,
274  new GeoField
275  (
276  IOobject
277  (
278  name,
279  mesh.thisDb().time().timeName(),
280  mesh.thisDb(),
281  IOobject::NO_READ,
282  IOobject::AUTO_WRITE
283  ),
284  mesh,
285  fieldDict
286  )
287  );
288 
290  // fields[i].write();
291  }
292  }
293  else
294  {
295  // Have mesh so just try to load
296  forAll(masterNames, i)
297  {
298  const word& name = masterNames[i];
299  IOobject& io = *objects[name];
300  io.writeOpt() = IOobject::AUTO_WRITE;
301 
302  // Load field
303  fields.set(i, new GeoField(io, mesh));
304  }
305  }
306 }
307 
308 
309 int main(int argc, char *argv[])
310 {
311  #include "addRegionOption.H"
312  #include "addOverwriteOption.H"
313 
314  // Include explicit constant options, have zero from time range
315  timeSelector::addOptions();
316 
317  #include "setRootCase.H"
318 
319  if (env("FOAM_SIGFPE"))
320  {
322  << "Detected floating point exception trapping (FOAM_SIGFPE)."
323  << " This might give" << nl
324  << " problems when mapping fields. Switch it off in case"
325  << " of problems." << endl;
326  }
327 
328 
329  // Create processor directory if non-existing
330  if (!Pstream::master() && !isDir(args.path()))
331  {
332  Pout<< "Creating case directory " << args.path() << endl;
333  mkDir(args.path());
334  }
335 
336 
337  // Make sure we do not use the master-only reading.
338  regIOobject::fileModificationChecking = regIOobject::timeStamp;
339  #include "createTime.H"
340  // Allow override of time
341  instantList times = timeSelector::selectIfPresent(runTime, args);
342  runTime.setTime(times[0], 0);
343 
344  runTime.functionObjects().off();
345 
346  word regionName = polyMesh::defaultRegion;
347  fileName meshSubDir;
348 
349  if (args.optionReadIfPresent("region", regionName))
350  {
351  meshSubDir = regionName/polyMesh::meshSubDir;
352  }
353  else
354  {
355  meshSubDir = polyMesh::meshSubDir;
356  }
357  Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
358 
359  const bool overwrite = args.optionFound("overwrite");
360 
361 
362  // Get time instance directory. Since not all processors have meshes
363  // just use the master one everywhere.
364 
365  fileName masterInstDir;
366  if (Pstream::master())
367  {
368  masterInstDir = runTime.findInstance(meshSubDir, "points");
369  }
370  Pstream::scatter(masterInstDir);
371 
372  // Check who has a mesh
373  const fileName meshPath = runTime.path()/masterInstDir/meshSubDir;
374 
375  Info<< "Found points in " << meshPath << nl << endl;
376 
377 
378  boolList haveMesh(Pstream::nProcs(), false);
379  haveMesh[Pstream::myProcNo()] = isDir(meshPath);
380  Pstream::gatherList(haveMesh);
381  Pstream::scatterList(haveMesh);
382  Info<< "Per processor mesh availability : " << haveMesh << endl;
383  const bool allHaveMesh = (findIndex(haveMesh, false) == -1);
384 
385  autoPtr<fvMesh> meshPtr = loadOrCreateMesh
386  (
387  IOobject
388  (
389  regionName,
390  masterInstDir,
391  runTime,
393  )
394  );
395 
396  fvMesh& mesh = meshPtr();
397 
398  // Print some statistics
399  Info<< "Before distribution:" << endl;
400  printMeshData(mesh);
401 
402 
403  labelList finalDecomp;
404 
405  // Create decompositionMethod and new decomposition
406  {
407  autoPtr<decompositionMethod> distributor
408  (
409  decompositionMethod::NewDistributor
410  (
411  decompositionMethod::decomposeParDict(runTime)
412  )
413  );
414 
415  finalDecomp = distributor().decompose(mesh, mesh.cellCentres());
416  }
417 
418  // Dump decomposition to volScalarField
419  if (!overwrite)
420  {
421  writeDecomposition("decomposition", mesh, finalDecomp);
422  }
423 
424 
425  // Create 0 sized mesh to do all the generation of zero sized
426  // fields on processors that have zero sized meshes. Note that this is
427  // only necessary on master but since polyMesh construction with
428  // Pstream::parRun does parallel comms we have to do it on all
429  // processors
430  autoPtr<fvMeshSubset> subsetterPtr;
431 
432  if (!allHaveMesh)
433  {
434  // Find last non-processor patch.
435  const polyBoundaryMesh& patches = mesh.boundaryMesh();
436 
437  label nonProci = -1;
438 
439  forAll(patches, patchi)
440  {
441  if (isA<processorPolyPatch>(patches[patchi]))
442  {
443  break;
444  }
445  nonProci++;
446  }
447 
448  if (nonProci == -1)
449  {
451  << "Cannot find non-processor patch on processor "
452  << Pstream::myProcNo() << endl
453  << " Current patches:" << patches.names() << abort(FatalError);
454  }
455 
456  // Subset 0 cells, no parallel comms. This is used to create zero-sized
457  // fields.
458  subsetterPtr.reset(new fvMeshSubset(mesh));
459  subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProci, false);
460  }
461 
462 
463  // Get original objects (before incrementing time!)
464  IOobjectList objects(mesh, runTime.timeName());
465  // We don't want to map the decomposition (mapping already tested when
466  // mapping the cell centre field)
467  IOobjectList::iterator iter = objects.find("decomposition");
468  if (iter != objects.end())
469  {
470  objects.erase(iter);
471  }
472 
473 
474  // volFields
475 
476  PtrList<volScalarField> volScalarFields;
477  readFields
478  (
479  haveMesh,
480  mesh,
481  subsetterPtr,
482  objects,
483  volScalarFields
484  );
485 
486  PtrList<volVectorField> volVectorFields;
487  readFields
488  (
489  haveMesh,
490  mesh,
491  subsetterPtr,
492  objects,
493  volVectorFields
494  );
495 
496  PtrList<volSphericalTensorField> volSphereTensorFields;
497  readFields
498  (
499  haveMesh,
500  mesh,
501  subsetterPtr,
502  objects,
503  volSphereTensorFields
504  );
505 
506  PtrList<volSymmTensorField> volSymmTensorFields;
507  readFields
508  (
509  haveMesh,
510  mesh,
511  subsetterPtr,
512  objects,
513  volSymmTensorFields
514  );
515 
516  PtrList<volTensorField> volTensorFields;
517  readFields
518  (
519  haveMesh,
520  mesh,
521  subsetterPtr,
522  objects,
523  volTensorFields
524  );
525 
526 
527  // surfaceFields
528 
529  PtrList<surfaceScalarField> surfScalarFields;
530  readFields
531  (
532  haveMesh,
533  mesh,
534  subsetterPtr,
535  objects,
536  surfScalarFields
537  );
538 
539  PtrList<surfaceVectorField> surfVectorFields;
540  readFields
541  (
542  haveMesh,
543  mesh,
544  subsetterPtr,
545  objects,
546  surfVectorFields
547  );
548 
549  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
550  readFields
551  (
552  haveMesh,
553  mesh,
554  subsetterPtr,
555  objects,
556  surfSphereTensorFields
557  );
558 
559  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
560  readFields
561  (
562  haveMesh,
563  mesh,
564  subsetterPtr,
565  objects,
566  surfSymmTensorFields
567  );
568 
569  PtrList<surfaceTensorField> surfTensorFields;
570  readFields
571  (
572  haveMesh,
573  mesh,
574  subsetterPtr,
575  objects,
576  surfTensorFields
577  );
578 
579 
580  // pointFields
581 
582  pointMesh& pMesh =
583  const_cast<pointMesh&>(pointMesh::New(mesh));
584  PtrList<pointScalarField> pointScalarFields;
585  readFields
586  (
587  haveMesh,
588  pMesh,
589  subsetterPtr,
590  objects,
591  pointScalarFields
592  );
593 
594  PtrList<pointVectorField> pointVectorFields;
595  readFields
596  (
597  haveMesh,
598  pMesh,
599  subsetterPtr,
600  objects,
601  pointVectorFields
602  );
603 
604  PtrList<pointSphericalTensorField> pointSphereTensorFields;
605  readFields
606  (
607  haveMesh,
608  pMesh,
609  subsetterPtr,
610  objects,
611  pointSphereTensorFields
612  );
613 
614  PtrList<pointSymmTensorField> pointSymmTensorFields;
615  readFields
616  (
617  haveMesh,
618  pMesh,
619  subsetterPtr,
620  objects,
621  pointSymmTensorFields
622  );
623 
624  PtrList<pointTensorField> pointTensorFields;
625  readFields
626  (
627  haveMesh,
628  pMesh,
629  subsetterPtr,
630  objects,
631  pointTensorFields
632  );
633 
634  // Debugging: Create additional volField that will be mapped.
635  // Used to test correctness of mapping
636  // volVectorField mapCc("mapCc", 1*mesh.C());
637 
638  // Mesh distribution engine
639  fvMeshDistribute distributor(mesh);
640 
641  // Pout<< "Wanted distribution:"
642  // << distributor.countCells(finalDecomp) << nl << endl;
643 
644  // Do actual sending/receiving of mesh
645  autoPtr<polyDistributionMap> map = distributor.distribute(finalDecomp);
646 
648  // map().distributeFaceData(faceCc);
649 
650 
651  // Print some statistics
652  Info<< "After distribution:" << endl;
653  printMeshData(mesh);
654 
655 
656  if (!overwrite)
657  {
658  runTime++;
659  }
660  else
661  {
662  mesh.setInstance(masterInstDir);
663  }
664  Info<< "Writing redistributed mesh to " << runTime.timeName() << nl << endl;
665  mesh.write();
666 
667 
668  // Print nice message
669  // ~~~~~~~~~~~~~~~~~~
670 
671  // Determine which processors remain so we can print nice final message.
672  labelList nFaces(Pstream::nProcs());
673  nFaces[Pstream::myProcNo()] = mesh.nFaces();
674  Pstream::gatherList(nFaces);
675  Pstream::scatterList(nFaces);
676 
677  Info<< nl
678  << "You can pick up the redecomposed mesh from the polyMesh directory"
679  << " in " << runTime.timeName() << "." << nl
680  << "If you redecomposed the mesh to less processors you can delete"
681  << nl
682  << "the processor directories with 0 sized meshes in them." << nl
683  << "Below is a sample set of commands to do this."
684  << " Take care when issuing these" << nl
685  << "commands." << nl << endl;
686 
687  forAll(nFaces, proci)
688  {
689  fileName procDir = "processor" + name(proci);
690 
691  if (nFaces[proci] == 0)
692  {
693  Info<< " rm -r " << procDir.c_str() << nl;
694  }
695  else
696  {
697  fileName timeDir = procDir/runTime.timeName()/meshSubDir;
698  fileName constDir = procDir/runTime.constant()/meshSubDir;
699 
700  Info<< " rm -r " << constDir.c_str() << nl
701  << " mv " << timeDir.c_str()
702  << ' ' << constDir.c_str() << nl;
703  }
704  }
705  Info<< endl;
706 
707 
708  Info<< "End\n" << endl;
709 
710  return 0;
711 }
712 
713 
714 // ************************************************************************* //
List< instant > instantList
List of instants.
Definition: instantList.H:42
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Inter-processor communication reduction functions.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
static fvMesh * meshPtr
Definition: globalFoam.H:52
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
const dimensionSet dimless
fvMesh & mesh
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
mkDir(pdfPath)
Load or create (0 size) a mesh. Used in distributing meshes to a larger number of processors...
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const char nl
Definition: Ostream.H:260
objects
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName path() const
Return the path to the caseName.
Definition: argListI.H:66
List< word > wordList
A List of words.
Definition: fileName.H:54
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
messageStream Info
bool env(const word &)
Return true if environment variable of given name is defined.
Definition: POSIX.C:91
Foam::argList args(argc, argv)