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-2021 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 "mapDistributePolyMesh.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 
404  IOdictionary decompositionDict
405  (
406  IOobject
407  (
408  "decomposeParDict",
409  runTime.system(),
410  mesh,
411  IOobject::MUST_READ_IF_MODIFIED,
412  IOobject::NO_WRITE
413  )
414  );
415 
416  labelList finalDecomp;
417 
418  // Create decompositionMethod and new decomposition
419  {
420  autoPtr<decompositionMethod> decomposer
421  (
423  (
424  decompositionDict
425  )
426  );
427 
428  if (!decomposer().parallelAware())
429  {
431  << "You have selected decomposition method "
432  << decomposer().typeName
433  << " which does" << endl
434  << "not synchronise the decomposition across"
435  << " processor patches." << endl
436  << " You might want to select a decomposition method which"
437  << " is aware of this. Continuing."
438  << endl;
439  }
440 
441  finalDecomp = decomposer().decompose(mesh, mesh.cellCentres());
442  }
443 
444  // Dump decomposition to volScalarField
445  if (!overwrite)
446  {
447  writeDecomposition("decomposition", mesh, finalDecomp);
448  }
449 
450 
451  // Create 0 sized mesh to do all the generation of zero sized
452  // fields on processors that have zero sized meshes. Note that this is
453  // only necessary on master but since polyMesh construction with
454  // Pstream::parRun does parallel comms we have to do it on all
455  // processors
456  autoPtr<fvMeshSubset> subsetterPtr;
457 
458  if (!allHaveMesh)
459  {
460  // Find last non-processor patch.
461  const polyBoundaryMesh& patches = mesh.boundaryMesh();
462 
463  label nonProci = -1;
464 
465  forAll(patches, patchi)
466  {
467  if (isA<processorPolyPatch>(patches[patchi]))
468  {
469  break;
470  }
471  nonProci++;
472  }
473 
474  if (nonProci == -1)
475  {
477  << "Cannot find non-processor patch on processor "
478  << Pstream::myProcNo() << endl
479  << " Current patches:" << patches.names() << abort(FatalError);
480  }
481 
482  // Subset 0 cells, no parallel comms. This is used to create zero-sized
483  // fields.
484  subsetterPtr.reset(new fvMeshSubset(mesh));
485  subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProci, false);
486  }
487 
488 
489  // Get original objects (before incrementing time!)
490  IOobjectList objects(mesh, runTime.timeName());
491  // We don't want to map the decomposition (mapping already tested when
492  // mapping the cell centre field)
493  IOobjectList::iterator iter = objects.find("decomposition");
494  if (iter != objects.end())
495  {
496  objects.erase(iter);
497  }
498 
499 
500  // volFields
501 
502  PtrList<volScalarField> volScalarFields;
503  readFields
504  (
505  haveMesh,
506  mesh,
507  subsetterPtr,
508  objects,
509  volScalarFields
510  );
511 
512  PtrList<volVectorField> volVectorFields;
513  readFields
514  (
515  haveMesh,
516  mesh,
517  subsetterPtr,
518  objects,
519  volVectorFields
520  );
521 
522  PtrList<volSphericalTensorField> volSphereTensorFields;
523  readFields
524  (
525  haveMesh,
526  mesh,
527  subsetterPtr,
528  objects,
529  volSphereTensorFields
530  );
531 
532  PtrList<volSymmTensorField> volSymmTensorFields;
533  readFields
534  (
535  haveMesh,
536  mesh,
537  subsetterPtr,
538  objects,
539  volSymmTensorFields
540  );
541 
542  PtrList<volTensorField> volTensorFields;
543  readFields
544  (
545  haveMesh,
546  mesh,
547  subsetterPtr,
548  objects,
549  volTensorFields
550  );
551 
552 
553  // surfaceFields
554 
555  PtrList<surfaceScalarField> surfScalarFields;
556  readFields
557  (
558  haveMesh,
559  mesh,
560  subsetterPtr,
561  objects,
562  surfScalarFields
563  );
564 
565  PtrList<surfaceVectorField> surfVectorFields;
566  readFields
567  (
568  haveMesh,
569  mesh,
570  subsetterPtr,
571  objects,
572  surfVectorFields
573  );
574 
575  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
576  readFields
577  (
578  haveMesh,
579  mesh,
580  subsetterPtr,
581  objects,
582  surfSphereTensorFields
583  );
584 
585  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
586  readFields
587  (
588  haveMesh,
589  mesh,
590  subsetterPtr,
591  objects,
592  surfSymmTensorFields
593  );
594 
595  PtrList<surfaceTensorField> surfTensorFields;
596  readFields
597  (
598  haveMesh,
599  mesh,
600  subsetterPtr,
601  objects,
602  surfTensorFields
603  );
604 
605 
606  // pointFields
607 
608  pointMesh& pMesh =
609  const_cast<pointMesh&>(pointMesh::New(mesh));
610  PtrList<pointScalarField> pointScalarFields;
611  readFields
612  (
613  haveMesh,
614  pMesh,
615  subsetterPtr,
616  objects,
617  pointScalarFields
618  );
619 
620  PtrList<pointVectorField> pointVectorFields;
621  readFields
622  (
623  haveMesh,
624  pMesh,
625  subsetterPtr,
626  objects,
627  pointVectorFields
628  );
629 
630  PtrList<pointSphericalTensorField> pointSphereTensorFields;
631  readFields
632  (
633  haveMesh,
634  pMesh,
635  subsetterPtr,
636  objects,
637  pointSphereTensorFields
638  );
639 
640  PtrList<pointSymmTensorField> pointSymmTensorFields;
641  readFields
642  (
643  haveMesh,
644  pMesh,
645  subsetterPtr,
646  objects,
647  pointSymmTensorFields
648  );
649 
650  PtrList<pointTensorField> pointTensorFields;
651  readFields
652  (
653  haveMesh,
654  pMesh,
655  subsetterPtr,
656  objects,
657  pointTensorFields
658  );
659 
660  // Debugging: Create additional volField that will be mapped.
661  // Used to test correctness of mapping
662  // volVectorField mapCc("mapCc", 1*mesh.C());
663 
664  // Mesh distribution engine
665  fvMeshDistribute distributor(mesh);
666 
667  // Pout<< "Wanted distribution:"
668  // << distributor.countCells(finalDecomp) << nl << endl;
669 
670  // Do actual sending/receiving of mesh
671  autoPtr<mapDistributePolyMesh> map = distributor.distribute(finalDecomp);
672 
674  // map().distributeFaceData(faceCc);
675 
676 
677  // Print some statistics
678  Info<< "After distribution:" << endl;
679  printMeshData(mesh);
680 
681 
682  if (!overwrite)
683  {
684  runTime++;
685  }
686  else
687  {
688  mesh.setInstance(masterInstDir);
689  }
690  Info<< "Writing redistributed mesh to " << runTime.timeName() << nl << endl;
691  mesh.write();
692 
693 
694  // Print nice message
695  // ~~~~~~~~~~~~~~~~~~
696 
697  // Determine which processors remain so we can print nice final message.
698  labelList nFaces(Pstream::nProcs());
699  nFaces[Pstream::myProcNo()] = mesh.nFaces();
700  Pstream::gatherList(nFaces);
701  Pstream::scatterList(nFaces);
702 
703  Info<< nl
704  << "You can pick up the redecomposed mesh from the polyMesh directory"
705  << " in " << runTime.timeName() << "." << nl
706  << "If you redecomposed the mesh to less processors you can delete"
707  << nl
708  << "the processor directories with 0 sized meshes in them." << nl
709  << "Below is a sample set of commands to do this."
710  << " Take care when issuing these" << nl
711  << "commands." << nl << endl;
712 
713  forAll(nFaces, proci)
714  {
715  fileName procDir = "processor" + name(proci);
716 
717  if (nFaces[proci] == 0)
718  {
719  Info<< " rm -r " << procDir.c_str() << nl;
720  }
721  else
722  {
723  fileName timeDir = procDir/runTime.timeName()/meshSubDir;
724  fileName constDir = procDir/runTime.constant()/meshSubDir;
725 
726  Info<< " rm -r " << constDir.c_str() << nl
727  << " mv " << timeDir.c_str()
728  << ' ' << constDir.c_str() << nl;
729  }
730  }
731  Info<< endl;
732 
733 
734  Info<< "End\n" << endl;
735 
736  return 0;
737 }
738 
739 
740 // ************************************************************************* //
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
virtual Ostream & write(const char)=0
Write character.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
static fvMesh * meshPtr
Definition: globalFoam.H:52
engineTime & runTime
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
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:57
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
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,.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
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)