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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 Usage
35 
36  - redistributePar [OPTION]
37 
38  \param -region regionName \n
39  Distribute named region.
40 
41  \param -decompose \n
42  Remove any existing \a processor subdirectories and decomposes the
43  geometry. Equivalent to running without processor subdirectories.
44 
45  \param -reconstruct \n
46  Reconstruct geometry (like reconstructParMesh). Equivalent to setting
47  numberOfSubdomains 1 in the decomposeParDict.
48 
49  \param -newTimes \n
50  (in combination with -reconstruct) reconstruct only new times.
51 
52  \param -cellDist \n
53  Write the cell distribution as a labelList, for use with 'manual'
54  decomposition method or as a volScalarField for post-processing.
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #include "argList.H"
59 #include "Time.H"
60 #include "fvMesh.H"
61 #include "fvMeshTools.H"
62 #include "fvMeshDistribute.H"
63 #include "decompositionMethod.H"
64 #include "timeSelector.H"
65 #include "PstreamReduceOps.H"
66 #include "volFields.H"
67 #include "surfaceFields.H"
69 #include "IOobjectList.H"
70 #include "globalIndex.H"
71 #include "loadOrCreateMesh.H"
72 #include "processorFvPatchField.H"
74 #include "decompositionModel.H"
75 
79 #include "hexRef8Data.H"
80 
81 using namespace Foam;
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 // Tolerance (as fraction of the bounding box). Needs to be fairly lax since
86 // usually meshes get written with limited precision (6 digits)
87 static const scalar defaultMergeTol = 1e-6;
88 
89 
90 // Get merging distance when matching face centres
91 scalar getMergeDistance
92 (
93  const argList& args,
94  const Time& runTime,
95  const boundBox& bb
96 )
97 {
98  scalar mergeTol = defaultMergeTol;
99  args.optionReadIfPresent("mergeTol", mergeTol);
100 
101  scalar writeTol =
102  Foam::pow(scalar(10.0), -scalar(IOstream::defaultPrecision()));
103 
104  Info<< "Merge tolerance : " << mergeTol << nl
105  << "Write tolerance : " << writeTol << endl;
106 
107  if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
108  {
110  << "Your current settings specify ASCII writing with "
111  << IOstream::defaultPrecision() << " digits precision." << endl
112  << "Your merging tolerance (" << mergeTol << ") is finer than this."
113  << endl
114  << "Please change your writeFormat to binary"
115  << " or increase the writePrecision" << endl
116  << "or adjust the merge tolerance (-mergeTol)."
117  << exit(FatalError);
118  }
119 
120  scalar mergeDist = mergeTol * bb.mag();
121 
122  Info<< "Overall meshes bounding box : " << bb << nl
123  << "Relative tolerance : " << mergeTol << nl
124  << "Absolute matching distance : " << mergeDist << nl
125  << endl;
126 
127  return mergeDist;
128 }
129 
130 
131 void printMeshData(const polyMesh& mesh)
132 {
133  // Collect all data on master
134 
135  globalIndex globalCells(mesh.nCells());
136  labelListList patchNeiProcNo(Pstream::nProcs());
137  labelListList patchSize(Pstream::nProcs());
138  const labelList& pPatches = mesh.globalData().processorPatches();
139  patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
140  patchSize[Pstream::myProcNo()].setSize(pPatches.size());
141  forAll(pPatches, i)
142  {
143  const processorPolyPatch& ppp = refCast<const processorPolyPatch>
144  (
145  mesh.boundaryMesh()[pPatches[i]]
146  );
147  patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
148  patchSize[Pstream::myProcNo()][i] = ppp.size();
149  }
150  Pstream::gatherList(patchNeiProcNo);
151  Pstream::gatherList(patchSize);
152 
153 
154  // Print stats
155 
156  globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());
157 
158  label maxProcCells = 0;
159  label totProcFaces = 0;
160  label maxProcPatches = 0;
161  label totProcPatches = 0;
162  label maxProcFaces = 0;
163 
164  for (label procI = 0; procI < Pstream::nProcs(); procI++)
165  {
166  Info<< endl
167  << "Processor " << procI << nl
168  << " Number of cells = " << globalCells.localSize(procI)
169  << endl;
170 
171  label nProcFaces = 0;
172 
173  const labelList& nei = patchNeiProcNo[procI];
174 
175  forAll(patchNeiProcNo[procI], i)
176  {
177  Info<< " Number of faces shared with processor "
178  << patchNeiProcNo[procI][i] << " = " << patchSize[procI][i]
179  << endl;
180 
181  nProcFaces += patchSize[procI][i];
182  }
183 
184  Info<< " Number of processor patches = " << nei.size() << nl
185  << " Number of processor faces = " << nProcFaces << nl
186  << " Number of boundary faces = "
187  << globalBoundaryFaces.localSize(procI)-nProcFaces << endl;
188 
189  maxProcCells = max(maxProcCells, globalCells.localSize(procI));
190  totProcFaces += nProcFaces;
191  totProcPatches += nei.size();
192  maxProcPatches = max(maxProcPatches, nei.size());
193  maxProcFaces = max(maxProcFaces, nProcFaces);
194  }
195 
196  // Stats
197 
198  scalar avgProcCells = scalar(globalCells.size())/Pstream::nProcs();
199  scalar avgProcPatches = scalar(totProcPatches)/Pstream::nProcs();
200  scalar avgProcFaces = scalar(totProcFaces)/Pstream::nProcs();
201 
202  // In case of all faces on one processor. Just to avoid division by 0.
203  if (totProcPatches == 0)
204  {
205  avgProcPatches = 1;
206  }
207  if (totProcFaces == 0)
208  {
209  avgProcFaces = 1;
210  }
211 
212  Info<< nl
213  << "Number of processor faces = " << totProcFaces/2 << nl
214  << "Max number of cells = " << maxProcCells
215  << " (" << 100.0*(maxProcCells-avgProcCells)/avgProcCells
216  << "% above average " << avgProcCells << ")" << nl
217  << "Max number of processor patches = " << maxProcPatches
218  << " (" << 100.0*(maxProcPatches-avgProcPatches)/avgProcPatches
219  << "% above average " << avgProcPatches << ")" << nl
220  << "Max number of faces between processors = " << maxProcFaces
221  << " (" << 100.0*(maxProcFaces-avgProcFaces)/avgProcFaces
222  << "% above average " << avgProcFaces << ")" << nl
223  << endl;
224 }
225 
226 
227 // Debugging: write volScalarField with decomposition for post processing.
229 (
230  const word& name,
231  const fvMesh& mesh,
232  const labelList& decomp
233 )
234 {
235  // Write the decomposition as labelList for use with 'manual'
236  // decomposition method.
237  labelIOList cellDecomposition
238  (
239  IOobject
240  (
241  "cellDecomposition",
242  mesh.facesInstance(), // mesh read from facesInstance
243  mesh,
246  false
247  ),
248  decomp
249  );
250  cellDecomposition.write();
251 
252  Info<< "Writing wanted cell distribution to volScalarField " << name
253  << " for postprocessing purposes." << nl << endl;
254 
255  volScalarField procCells
256  (
257  IOobject
258  (
259  name,
260  mesh.time().timeName(),
261  mesh,
264  false // do not register
265  ),
266  mesh,
268  zeroGradientFvPatchScalarField::typeName
269  );
270 
271  forAll(procCells, cI)
272  {
273  procCells[cI] = decomp[cI];
274  }
275  procCells.write();
276 }
277 
278 
280 (
281  const Time& baseRunTime,
282  const fileName& decompDictFile, // optional location for decomposeParDict
283  const bool decompose, // decompose, i.e. read from undecomposed case
284  const fileName& proc0CaseName,
285  const fvMesh& mesh,
286  const bool writeCellDist,
287 
288  label& nDestProcs,
289  labelList& decomp
290 )
291 {
292  // Read decomposeParDict (on all processors)
294  (
295  mesh,
296  decompDictFile
297  );
298 
299  decompositionMethod& decomposer = method.decomposer();
300 
301  if (!decomposer.parallelAware())
302  {
304  << "You have selected decomposition method "
305  << decomposer.typeName
306  << " which does" << endl
307  << "not synchronise the decomposition across"
308  << " processor patches." << endl
309  << " You might want to select a decomposition method"
310  << " which is aware of this. Continuing."
311  << endl;
312  }
313 
314  if (Pstream::master() && decompose)
315  {
316  Info<< "Setting caseName to " << baseRunTime.caseName()
317  << " to read decomposeParDict" << endl;
318  const_cast<Time&>(mesh.time()).TimePaths::caseName() =
319  baseRunTime.caseName();
320  }
321 
322  scalarField cellWeights;
323  if (method.found("weightField"))
324  {
325  word weightName = method.lookup("weightField");
326 
327  volScalarField weights
328  (
329  IOobject
330  (
331  weightName,
332  mesh.time().timeName(),
333  mesh,
336  ),
337  mesh
338  );
339  cellWeights = weights.internalField();
340  }
341 
342  nDestProcs = decomposer.nDomains();
343  decomp = decomposer.decompose(mesh, cellWeights);
344 
345  if (Pstream::master() && decompose)
346  {
347  Info<< "Restoring caseName to " << proc0CaseName << endl;
348  const_cast<Time&>(mesh.time()).TimePaths::caseName() =
349  proc0CaseName;
350  }
351 
352  // Dump decomposition to volScalarField
353  if (writeCellDist)
354  {
355  // Note: on master make sure to write to processor0
356  if (decompose)
357  {
358  if (Pstream::master())
359  {
360  Info<< "Setting caseName to " << baseRunTime.caseName()
361  << " to write undecomposed cellDist" << endl;
362 
363  Time& tm = const_cast<Time&>(mesh.time());
364 
365  tm.TimePaths::caseName() = baseRunTime.caseName();
366  writeDecomposition("cellDist", mesh, decomp);
367  Info<< "Restoring caseName to " << proc0CaseName << endl;
368  tm.TimePaths::caseName() = proc0CaseName;
369  }
370  }
371  else
372  {
373  writeDecomposition("cellDist", mesh, decomp);
374  }
375  }
376 }
377 
378 
379 // Write addressing if decomposing (1 to many) or reconstructing (many to 1)
381 (
382  const bool decompose,
383  const fileName& meshSubDir,
384  const fvMesh& mesh,
385  const mapDistributePolyMesh& map
386 )
387 {
388  Info<< "Writing procAddressing files to " << mesh.facesInstance()
389  << endl;
390 
391  labelIOList cellMap
392  (
393  IOobject
394  (
395  "cellProcAddressing",
397  meshSubDir,
398  mesh,
400  ),
401  0
402  );
403 
405  (
406  IOobject
407  (
408  "faceProcAddressing",
410  meshSubDir,
411  mesh,
413  ),
414  0
415  );
416 
417  labelIOList pointMap
418  (
419  IOobject
420  (
421  "pointProcAddressing",
423  meshSubDir,
424  mesh,
426  ),
427  0
428  );
429 
430  labelIOList patchMap
431  (
432  IOobject
433  (
434  "boundaryProcAddressing",
436  meshSubDir,
437  mesh,
439  ),
440  0
441  );
442 
443  // Decomposing: see how cells moved from undecomposed case
444  if (decompose)
445  {
446  cellMap = identity(map.nOldCells());
447  map.distributeCellData(cellMap);
448 
449  faceMap = identity(map.nOldFaces());
450  {
451  const mapDistribute& faceDistMap = map.faceMap();
452 
453  if (faceDistMap.subHasFlip() || faceDistMap.constructHasFlip())
454  {
455  // Offset by 1
456  faceMap = faceMap + 1;
457  }
458  // Apply face flips
460  (
462  List<labelPair>(),
463  faceDistMap.constructSize(),
464  faceDistMap.subMap(),
465  faceDistMap.subHasFlip(),
466  faceDistMap.constructMap(),
467  faceDistMap.constructHasFlip(),
468  faceMap,
469  flipLabelOp()
470  );
471  }
472 
473  pointMap = identity(map.nOldPoints());
474  map.distributePointData(pointMap);
475 
476  patchMap = identity(map.oldPatchSizes().size());
477  const mapDistribute& patchDistMap = map.patchMap();
478  // Use explicit distribute since we need to provide a null value
479  // (for new patches) and this is the only call that allow us to
480  // provide one ...
482  (
484  List<labelPair>(),
485  patchDistMap.constructSize(),
486  patchDistMap.subMap(),
487  patchDistMap.subHasFlip(),
488  patchDistMap.constructMap(),
489  patchDistMap.constructHasFlip(),
490  patchMap,
491  eqOp<label>(),
492  flipOp(),
493  label(-1),
495  );
496  }
497  else // if (nDestProcs == 1)
498  {
499  cellMap = identity(mesh.nCells());
500  map.cellMap().reverseDistribute(map.nOldCells(), cellMap);
501 
503  {
504  const mapDistribute& faceDistMap = map.faceMap();
505 
506  if (faceDistMap.subHasFlip() || faceDistMap.constructHasFlip())
507  {
508  // Offset by 1
509  faceMap = faceMap + 1;
510  }
511 
513  (
515  List<labelPair>(),
516  map.nOldFaces(),
517  faceDistMap.constructMap(),
518  faceDistMap.constructHasFlip(),
519  faceDistMap.subMap(),
520  faceDistMap.subHasFlip(),
521  faceMap,
522  flipLabelOp()
523  );
524  }
525 
526  pointMap = identity(mesh.nPoints());
527  map.pointMap().reverseDistribute(map.nOldPoints(), pointMap);
528 
529  const mapDistribute& patchDistMap = map.patchMap();
530  patchMap = identity(mesh.boundaryMesh().size());
531  patchDistMap.reverseDistribute
532  (
533  map.oldPatchSizes().size(),
534  label(-1),
535  patchMap
536  );
537  }
538 
539  bool cellOk = cellMap.write();
540  bool faceOk = faceMap.write();
541  bool pointOk = pointMap.write();
542  bool patchOk = patchMap.write();
543 
544  if (!cellOk || !faceOk || !pointOk || !patchOk)
545  {
547  << "Failed to write " << cellMap.objectPath()
548  << ", " << faceMap.objectPath()
549  << ", " << pointMap.objectPath()
550  << ", " << patchMap.objectPath()
551  << endl;
552  }
553 }
554 
555 
556 
557 // Generic mesh-based field reading
558 template<class GeoField>
559 void readField
560 (
561  const IOobject& io,
562  const fvMesh& mesh,
563  const label i,
565 )
566 {
567  fields.set(i, new GeoField(io, mesh));
568 }
569 
570 
571 // Definition of readField for GeometricFields only
572 template<class Type, template<class> class PatchField, class GeoMesh>
573 void readField
574 (
575  const IOobject& io,
576  const fvMesh& mesh,
577  const label i,
579 )
580 {
581  fields.set
582  (
583  i,
585  );
586 }
587 
588 
589 // Read vol or surface fields
590 template<class GeoField>
591 void readFields
592 (
593  const boolList& haveMesh,
594  const fvMesh& mesh,
595  const autoPtr<fvMeshSubset>& subsetterPtr,
596  IOobjectList& allObjects,
598 )
599 {
600  //typedef GeometricField<T, fvPatchField, Mesh> fldType;
601 
602  // Get my objects of type
603  IOobjectList objects(allObjects.lookupClass(GeoField::typeName));
604  // Check that we all have all objects
605  wordList objectNames = objects.sortedNames();
606  // Get master names
607  wordList masterNames(objectNames);
608  Pstream::scatter(masterNames);
609 
610  if (haveMesh[Pstream::myProcNo()] && objectNames != masterNames)
611  {
613  << "differing fields of type " << GeoField::typeName
614  << " on processors." << endl
615  << "Master has:" << masterNames << endl
616  << Pstream::myProcNo() << " has:" << objectNames
617  << abort(FatalError);
618  }
619 
620  fields.setSize(masterNames.size());
621 
622  // Have master send all fields to processors that don't have a mesh
623  if (Pstream::master())
624  {
625  forAll(masterNames, i)
626  {
627  const word& name = masterNames[i];
628  IOobject& io = *objects[name];
630 
631  // Load field (but not oldTime)
632  readField(io, mesh, i, fields);
633  // Create zero sized field and send
634  if (subsetterPtr.valid())
635  {
636  tmp<GeoField> tsubfld = subsetterPtr().interpolate(fields[i]);
637 
638  // Send to all processors that don't have a mesh
639  for (label procI = 1; procI < Pstream::nProcs(); procI++)
640  {
641  if (!haveMesh[procI])
642  {
643  OPstream toProc(Pstream::blocking, procI);
644  toProc<< tsubfld();
645  }
646  }
647  }
648  }
649  }
650  else if (!haveMesh[Pstream::myProcNo()])
651  {
652  // Don't have mesh (nor fields). Receive empty field from master.
653 
654  forAll(masterNames, i)
655  {
656  const word& name = masterNames[i];
657 
658  // Receive field
660  dictionary fieldDict(fromMaster);
661 
662  fields.set
663  (
664  i,
665  new GeoField
666  (
667  IOobject
668  (
669  name,
670  mesh.time().timeName(),
671  mesh,
674  ),
675  mesh,
676  fieldDict
677  )
678  );
679 
681  //fields[i].write();
682  }
683  }
684  else
685  {
686  // Have mesh so just try to load
687  forAll(masterNames, i)
688  {
689  const word& name = masterNames[i];
690  IOobject& io = *objects[name];
692 
693  // Load field (but not oldtime)
694  readField(io, mesh, i, fields);
695  }
696  }
697 }
698 
699 
700 // Variant of GeometricField::correctBoundaryConditions that only
701 // evaluates selected patch fields
702 template<class GeoField, class CoupledPatchType>
704 {
706  (
707  mesh.objectRegistry::lookupClass<GeoField>()
708  );
709 
710  forAllIter(typename HashTable<GeoField*>, flds, iter)
711  {
712  GeoField& fld = *iter();
713 
714  typename GeoField::GeometricBoundaryField& bfld =
715  fld.boundaryField();
716  if
717  (
720  )
721  {
722  label nReq = Pstream::nRequests();
723 
724  forAll(bfld, patchi)
725  {
726  typename GeoField::PatchFieldType& pfld = bfld[patchi];
727 
728  //if (pfld.coupled())
729  //if (isA<CoupledPatchType>(pfld))
730  if (pfld.patch().coupled())
731  {
732  pfld.initEvaluate(Pstream::defaultCommsType);
733  }
734  }
735 
736  // Block for any outstanding requests
737  if
738  (
741  )
742  {
743  Pstream::waitRequests(nReq);
744  }
745 
746  forAll(bfld, patchi)
747  {
748  typename GeoField::PatchFieldType& pfld = bfld[patchi];
749 
750  //if (pfld.coupled())
751  //if (isA<CoupledPatchType>(pfld))
752  if (pfld.patch().coupled())
753  {
754  pfld.evaluate(Pstream::defaultCommsType);
755  }
756  }
757  }
759  {
760  const lduSchedule& patchSchedule =
761  fld.mesh().globalData().patchSchedule();
762 
763  forAll(patchSchedule, patchEvali)
764  {
765  label patchi = patchSchedule[patchEvali].patch;
766  typename GeoField::PatchFieldType& pfld = bfld[patchi];
767 
768  //if (pfld.coupled())
769  //if (isA<CoupledPatchType>(pfld))
770  if (pfld.patch().coupled())
771  {
772  if (patchSchedule[patchEvali].init)
773  {
774  pfld.initEvaluate(Pstream::scheduled);
775  }
776  else
777  {
778  pfld.evaluate(Pstream::scheduled);
779  }
780  }
781  }
782  }
783  else
784  {
786  << "Unsuported communications type "
788  << exit(FatalError);
789  }
790  }
791 }
792 
793 
794 // Inplace redistribute mesh and any fields
796 (
797  const Time& baseRunTime,
798  const scalar tolDim,
799  const boolList& haveMesh,
800  const fileName& meshSubDir,
801  const bool doReadFields,
802  const bool decompose, // decompose, i.e. read from undecomposed case
803  const bool overwrite,
804  const fileName& proc0CaseName,
805  const label nDestProcs,
806  const labelList& decomp,
807  const fileName& masterInstDir,
808  fvMesh& mesh
809 )
810 {
811  Time& runTime = const_cast<Time&>(mesh.time());
812 
814  //Info<< "Before distribution:" << endl;
815  //printMeshData(mesh);
816 
817 
818  PtrList<volScalarField> volScalarFields;
819  PtrList<volVectorField> volVectorFields;
820  PtrList<volSphericalTensorField> volSphereTensorFields;
821  PtrList<volSymmTensorField> volSymmTensorFields;
822  PtrList<volTensorField> volTensorFields;
823 
824  PtrList<surfaceScalarField> surfScalarFields;
825  PtrList<surfaceVectorField> surfVectorFields;
826  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
827  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
828  PtrList<surfaceTensorField> surfTensorFields;
829 
833  PtrList<DimensionedField<symmTensor, volMesh> > dimSymmTensorFields;
835 
836 
837  if (doReadFields)
838  {
839  // Create 0 sized mesh to do all the generation of zero sized
840  // fields on processors that have zero sized meshes. Note that this is
841  // only nessecary on master but since polyMesh construction with
842  // Pstream::parRun does parallel comms we have to do it on all
843  // processors
844  autoPtr<fvMeshSubset> subsetterPtr;
845 
846  const bool allHaveMesh = (findIndex(haveMesh, false) == -1);
847  if (!allHaveMesh)
848  {
849  // Find last non-processor patch.
851 
852  label nonProcI = -1;
853 
854  forAll(patches, patchI)
855  {
856  if (isA<processorPolyPatch>(patches[patchI]))
857  {
858  break;
859  }
860  nonProcI++;
861  }
862 
863  if (nonProcI == -1)
864  {
866  << "Cannot find non-processor patch on processor "
867  << Pstream::myProcNo() << endl
868  << " Current patches:" << patches.names()
869  << abort(FatalError);
870  }
871 
872  // Subset 0 cells, no parallel comms. This is used to create
873  // zero-sized fields.
874  subsetterPtr.reset(new fvMeshSubset(mesh));
875  subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProcI, false);
876  }
877 
878 
879  // Get original objects (before incrementing time!)
880  if (Pstream::master() && decompose)
881  {
882  runTime.TimePaths::caseName() = baseRunTime.caseName();
883  }
884  IOobjectList objects(mesh, runTime.timeName());
885  if (Pstream::master() && decompose)
886  {
887  runTime.TimePaths::caseName() = proc0CaseName;
888  }
889 
890  Info<< "From time " << runTime.timeName()
891  << " have objects:" << objects.names() << endl;
892 
893  // We don't want to map the decomposition (mapping already tested when
894  // mapping the cell centre field)
895  IOobjectList::iterator iter = objects.find("cellDist");
896  if (iter != objects.end())
897  {
898  objects.erase(iter);
899  }
900 
901 
902  // volFields
903 
904  if (Pstream::master() && decompose)
905  {
906  runTime.TimePaths::caseName() = baseRunTime.caseName();
907  }
908  readFields
909  (
910  haveMesh,
911  mesh,
912  subsetterPtr,
913  objects,
914  volScalarFields
915  );
916 
917  readFields
918  (
919  haveMesh,
920  mesh,
921  subsetterPtr,
922  objects,
923  volVectorFields
924  );
925 
926  readFields
927  (
928  haveMesh,
929  mesh,
930  subsetterPtr,
931  objects,
932  volSphereTensorFields
933  );
934 
935  readFields
936  (
937  haveMesh,
938  mesh,
939  subsetterPtr,
940  objects,
941  volSymmTensorFields
942  );
943 
944  readFields
945  (
946  haveMesh,
947  mesh,
948  subsetterPtr,
949  objects,
950  volTensorFields
951  );
952 
953 
954  // surfaceFields
955 
956  readFields
957  (
958  haveMesh,
959  mesh,
960  subsetterPtr,
961  objects,
962  surfScalarFields
963  );
964 
965  readFields
966  (
967  haveMesh,
968  mesh,
969  subsetterPtr,
970  objects,
971  surfVectorFields
972  );
973 
974  readFields
975  (
976  haveMesh,
977  mesh,
978  subsetterPtr,
979  objects,
980  surfSphereTensorFields
981  );
982 
983  readFields
984  (
985  haveMesh,
986  mesh,
987  subsetterPtr,
988  objects,
989  surfSymmTensorFields
990  );
991 
992  readFields
993  (
994  haveMesh,
995  mesh,
996  subsetterPtr,
997  objects,
998  surfTensorFields
999  );
1000 
1001 
1002  // Dimensioned internal fields
1003  readFields
1004  (
1005  haveMesh,
1006  mesh,
1007  subsetterPtr,
1008  objects,
1009  dimScalarFields
1010  );
1011 
1012  readFields
1013  (
1014  haveMesh,
1015  mesh,
1016  subsetterPtr,
1017  objects,
1018  dimVectorFields
1019  );
1020 
1021  readFields
1022  (
1023  haveMesh,
1024  mesh,
1025  subsetterPtr,
1026  objects,
1027  dimSphereTensorFields
1028  );
1029 
1030  readFields
1031  (
1032  haveMesh,
1033  mesh,
1034  subsetterPtr,
1035  objects,
1036  dimSymmTensorFields
1037  );
1038 
1039  readFields
1040  (
1041  haveMesh,
1042  mesh,
1043  subsetterPtr,
1044  objects,
1045  dimTensorFields
1046  );
1047 
1048 
1049 
1050  if (Pstream::master() && decompose)
1051  {
1052  runTime.TimePaths::caseName() = proc0CaseName;
1053  }
1054  }
1055 
1056 
1057  // Mesh distribution engine
1058  fvMeshDistribute distributor(mesh, tolDim);
1059 
1060  // Do all the distribution of mesh and fields
1061  autoPtr<mapDistributePolyMesh> rawMap = distributor.distribute(decomp);
1062 
1063  // Print some statistics
1064  Info<< "After distribution:" << endl;
1066 
1067  // Get other side of processor boundaries
1069  <
1072  >(mesh);
1074  <
1077  >(mesh);
1079  <
1082  >(mesh);
1084  <
1087  >(mesh);
1089  <
1092  >(mesh);
1093  // No update surface fields
1094 
1095 
1096  // Set the minimum write precision
1098 
1099 
1100  if (!overwrite)
1101  {
1102  runTime++;
1103  mesh.setInstance(runTime.timeName());
1104  }
1105  else
1106  {
1107  mesh.setInstance(masterInstDir);
1108  }
1109 
1110 
1112  (
1113  IOobject
1114  (
1115  "procAddressing",
1116  mesh.facesInstance(),
1117  meshSubDir,
1118  mesh,
1121  )
1122  );
1123  map.transfer(rawMap());
1124 
1125 
1126  if (nDestProcs == 1)
1127  {
1128  if (Pstream::master())
1129  {
1130  Info<< "Setting caseName to " << baseRunTime.caseName()
1131  << " to write reconstructed mesh and fields." << endl;
1132  runTime.TimePaths::caseName() = baseRunTime.caseName();
1133 
1134  mesh.write();
1135 
1136  // Now we've written all. Reset caseName on master
1137  Info<< "Restoring caseName to " << proc0CaseName << endl;
1138  runTime.TimePaths::caseName() = proc0CaseName;
1139  }
1140  }
1141  else
1142  {
1143  mesh.write();
1144  }
1145  Info<< "Written redistributed mesh to " << mesh.facesInstance() << nl
1146  << endl;
1147 
1148 
1149  if (decompose || nDestProcs == 1)
1150  {
1151  writeProcAddressing(decompose, meshSubDir, mesh, map);
1152  }
1153 
1154 
1155  // Refinement data
1156  {
1157 
1158  // Read refinement data
1159  if (Pstream::master() && decompose)
1160  {
1161  runTime.TimePaths::caseName() = baseRunTime.caseName();
1162  }
1163  IOobject io
1164  (
1165  "dummy",
1166  mesh.facesInstance(),
1168  mesh,
1171  false
1172  );
1173 
1174  hexRef8Data refData(io);
1175  if (Pstream::master() && decompose)
1176  {
1177  runTime.TimePaths::caseName() = proc0CaseName;
1178  }
1179  // Make sure all processors have valid data (since only some will
1180  // read)
1181  refData.sync(io);
1182 
1183 
1184  // Distribute
1185  refData.distribute(map);
1186 
1187  if (nDestProcs == 1)
1188  {
1189  if (Pstream::master())
1190  {
1191  Info<< "Setting caseName to " << baseRunTime.caseName()
1192  << " to write reconstructed refinement data." << endl;
1193  runTime.TimePaths::caseName() = baseRunTime.caseName();
1194 
1195  refData.write();
1196 
1197  // Now we've written all. Reset caseName on master
1198  Info<< "Restoring caseName to " << proc0CaseName << endl;
1199  runTime.TimePaths::caseName() = proc0CaseName;
1200  }
1201  }
1202  else
1203  {
1204  refData.write();
1205  }
1206  }
1207 
1208 
1210  (
1211  new mapDistributePolyMesh(map.xfer())
1212  );
1213 }
1214 
1215 
1216 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1217 //
1218 // Field Mapping
1219 //
1220 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1221 
1223 (
1224  const autoPtr<fvMesh>& baseMeshPtr,
1225  const fvMesh& mesh,
1226  const labelList& cellProcAddressing,
1228  const labelList& pointProcAddressing,
1229  const labelList& boundaryProcAddressing
1230 )
1231 {
1232  // Send addressing to master
1233  labelListList cellAddressing(Pstream::nProcs());
1234  cellAddressing[Pstream::myProcNo()] = cellProcAddressing;
1235  Pstream::gatherList(cellAddressing);
1236 
1237  labelListList faceAddressing(Pstream::nProcs());
1238  faceAddressing[Pstream::myProcNo()] = faceProcAddressing;
1239  Pstream::gatherList(faceAddressing);
1240 
1241  labelListList pointAddressing(Pstream::nProcs());
1242  pointAddressing[Pstream::myProcNo()] = pointProcAddressing;
1243  Pstream::gatherList(pointAddressing);
1244 
1245  labelListList boundaryAddressing(Pstream::nProcs());
1246  {
1247  // Remove -1 entries
1248  DynamicList<label> patchProcAddressing(boundaryProcAddressing.size());
1249  forAll(boundaryProcAddressing, i)
1250  {
1251  if (boundaryProcAddressing[i] != -1)
1252  {
1253  patchProcAddressing.append(boundaryProcAddressing[i]);
1254  }
1255  }
1256  boundaryAddressing[Pstream::myProcNo()] = patchProcAddressing;
1257  Pstream::gatherList(boundaryAddressing);
1258  }
1259 
1260 
1262 
1263  if (baseMeshPtr.valid() && baseMeshPtr().nCells()) //baseMeshPtr.valid())
1264  {
1265  const fvMesh& baseMesh = baseMeshPtr();
1266 
1267  labelListList cellSubMap(Pstream::nProcs());
1268  cellSubMap[Pstream::masterNo()] = identity(mesh.nCells());
1269 
1270  mapDistribute cellMap
1271  (
1272  baseMesh.nCells(),
1273  cellSubMap.xfer(),
1274  cellAddressing.xfer()
1275  );
1276 
1277  labelListList faceSubMap(Pstream::nProcs());
1278  faceSubMap[Pstream::masterNo()] = identity(mesh.nFaces());
1279 
1281  (
1282  baseMesh.nFaces(),
1283  faceSubMap.xfer(),
1284  faceAddressing.xfer(),
1285  false, //subHasFlip
1286  true //constructHasFlip
1287  );
1288 
1289  labelListList pointSubMap(Pstream::nProcs());
1290  pointSubMap[Pstream::masterNo()] = identity(mesh.nPoints());
1291 
1292  mapDistribute pointMap
1293  (
1294  baseMesh.nPoints(),
1295  pointSubMap.xfer(),
1296  pointAddressing.xfer()
1297  );
1298 
1299  labelListList patchSubMap(Pstream::nProcs());
1300  // Send (filtered) patches to master
1301  patchSubMap[Pstream::masterNo()] =
1302  boundaryAddressing[Pstream::myProcNo()];
1303 
1304  mapDistribute patchMap
1305  (
1306  baseMesh.boundaryMesh().size(),
1307  patchSubMap.xfer(),
1308  boundaryAddressing.xfer()
1309  );
1310 
1311  const label nOldPoints = mesh.nPoints();
1312  const label nOldFaces = mesh.nFaces();
1313  const label nOldCells = mesh.nCells();
1314 
1315  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1316  labelList oldPatchStarts(pbm.size());
1317  labelList oldPatchNMeshPoints(pbm.size());
1318  forAll(pbm, patchI)
1319  {
1320  oldPatchStarts[patchI] = pbm[patchI].start();
1321  oldPatchNMeshPoints[patchI] = pbm[patchI].nPoints();
1322  }
1323 
1324  mapPtr.reset
1325  (
1327  (
1328  nOldPoints,
1329  nOldFaces,
1330  nOldCells,
1331  oldPatchStarts.xfer(),
1332  oldPatchNMeshPoints.xfer(),
1333  pointMap.xfer(),
1334  faceMap.xfer(),
1335  cellMap.xfer(),
1336  patchMap.xfer()
1337  )
1338  );
1339  }
1340  else
1341  {
1342  labelListList cellSubMap(Pstream::nProcs());
1343  cellSubMap[Pstream::masterNo()] = identity(mesh.nCells());
1344  labelListList cellConstructMap(Pstream::nProcs());
1345 
1346  mapDistribute cellMap
1347  (
1348  0,
1349  cellSubMap.xfer(),
1350  cellConstructMap.xfer()
1351  );
1352 
1353  labelListList faceSubMap(Pstream::nProcs());
1354  faceSubMap[Pstream::masterNo()] = identity(mesh.nFaces());
1355  labelListList faceConstructMap(Pstream::nProcs());
1356 
1358  (
1359  0,
1360  faceSubMap.xfer(),
1361  faceConstructMap.xfer(),
1362  false, //subHasFlip
1363  true //constructHasFlip
1364  );
1365 
1366  labelListList pointSubMap(Pstream::nProcs());
1367  pointSubMap[Pstream::masterNo()] = identity(mesh.nPoints());
1368  labelListList pointConstructMap(Pstream::nProcs());
1369 
1370  mapDistribute pointMap
1371  (
1372  0,
1373  pointSubMap.xfer(),
1374  pointConstructMap.xfer()
1375  );
1376 
1377  labelListList patchSubMap(Pstream::nProcs());
1378  // Send (filtered) patches to master
1379  patchSubMap[Pstream::masterNo()] =
1380  boundaryAddressing[Pstream::myProcNo()];
1381  labelListList patchConstructMap(Pstream::nProcs());
1382 
1383  mapDistribute patchMap
1384  (
1385  0,
1386  patchSubMap.xfer(),
1387  patchConstructMap.xfer()
1388  );
1389 
1390  const label nOldPoints = mesh.nPoints();
1391  const label nOldFaces = mesh.nFaces();
1392  const label nOldCells = mesh.nCells();
1393 
1394  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1395  labelList oldPatchStarts(pbm.size());
1396  labelList oldPatchNMeshPoints(pbm.size());
1397  forAll(pbm, patchI)
1398  {
1399  oldPatchStarts[patchI] = pbm[patchI].start();
1400  oldPatchNMeshPoints[patchI] = pbm[patchI].nPoints();
1401  }
1402 
1403  mapPtr.reset
1404  (
1406  (
1407  nOldPoints,
1408  nOldFaces,
1409  nOldCells,
1410  oldPatchStarts.xfer(),
1411  oldPatchNMeshPoints.xfer(),
1412  pointMap.xfer(),
1413  faceMap.xfer(),
1414  cellMap.xfer(),
1415  patchMap.xfer()
1416  )
1417  );
1418  }
1419 
1420  return mapPtr;
1421 }
1423 
1424 void readProcAddressing
1425 (
1426  const fvMesh& mesh,
1427  const autoPtr<fvMesh>& baseMeshPtr,
1429 )
1430 {
1431  //IOobject io
1432  //(
1433  // "procAddressing",
1434  // mesh.facesInstance(),
1435  // polyMesh::meshSubDir,
1436  // mesh,
1437  // IOobject::MUST_READ
1438  //);
1439  //if (io.headerOk())
1440  //{
1441  // Pout<< "Reading addressing from " << io.name() << " at "
1442  // << mesh.facesInstance() << nl << endl;
1443  // distMap.clear();
1444  // distMap.reset(new IOmapDistributePolyMesh(io));
1445  //}
1446  //else
1447  {
1448  Info<< "Reading addressing from procXXXAddressing at "
1449  << mesh.facesInstance() << nl << endl;
1450  labelIOList cellProcAddressing
1451  (
1452  IOobject
1453  (
1454  "cellProcAddressing",
1455  mesh.facesInstance(),
1457  mesh,
1459  ),
1460  labelList(0)
1461  );
1463  (
1464  IOobject
1465  (
1466  "faceProcAddressing",
1467  mesh.facesInstance(),
1469  mesh,
1471  ),
1472  labelList(0)
1473  );
1474  labelIOList pointProcAddressing
1475  (
1476  IOobject
1477  (
1478  "pointProcAddressing",
1479  mesh.facesInstance(),
1481  mesh,
1483  ),
1484  labelList(0)
1485  );
1486  labelIOList boundaryProcAddressing
1487  (
1488  IOobject
1489  (
1490  "boundaryProcAddressing",
1491  mesh.facesInstance(),
1493  mesh,
1495  ),
1496  labelList(0)
1497  );
1498 
1499 
1500  if
1501  (
1502  mesh.nCells() != cellProcAddressing.size()
1503  || mesh.nPoints() != pointProcAddressing.size()
1504  || mesh.nFaces() != faceProcAddressing.size()
1505  || mesh.boundaryMesh().size() != boundaryProcAddressing.size()
1506  )
1507  {
1509  << "Read addressing inconsistent with mesh sizes" << nl
1510  << "cells:" << mesh.nCells()
1511  << " addressing:" << cellProcAddressing.objectPath()
1512  << " size:" << cellProcAddressing.size() << nl
1513  << "faces:" << mesh.nFaces()
1514  << " addressing:" << faceProcAddressing.objectPath()
1515  << " size:" << faceProcAddressing.size() << nl
1516  << "points:" << mesh.nPoints()
1517  << " addressing:" << pointProcAddressing.objectPath()
1518  << " size:" << pointProcAddressing.size()
1519  << "patches:" << mesh.boundaryMesh().size()
1520  << " addressing:" << boundaryProcAddressing.objectPath()
1521  << " size:" << boundaryProcAddressing.size()
1522  << exit(FatalError);
1523  }
1524 
1525  distMap.clear();
1526  distMap = createReconstructMap
1527  (
1528  baseMeshPtr,
1529  mesh,
1530  cellProcAddressing,
1532  pointProcAddressing,
1533  boundaryProcAddressing
1534  );
1535  }
1536 }
1538 
1540 (
1541  const parFvFieldReconstructor& fvReconstructor,
1542  const IOobjectList& objects,
1543  const HashSet<word>& selectedFields
1544 )
1545 {
1546  // Dimensioned fields
1547 
1548  fvReconstructor.reconstructFvVolumeInternalFields<scalar>
1549  (
1550  objects,
1551  selectedFields
1552  );
1553  fvReconstructor.reconstructFvVolumeInternalFields<vector>
1554  (
1555  objects,
1556  selectedFields
1557  );
1559  (
1560  objects,
1561  selectedFields
1562  );
1564  (
1565  objects,
1566  selectedFields
1567  );
1568  fvReconstructor.reconstructFvVolumeInternalFields<tensor>
1569  (
1570  objects,
1571  selectedFields
1572  );
1573 
1574 
1575  // volFields
1576 
1577  fvReconstructor.reconstructFvVolumeFields<scalar>
1578  (
1579  objects,
1580  selectedFields
1581  );
1582  fvReconstructor.reconstructFvVolumeFields<vector>
1583  (
1584  objects,
1585  selectedFields
1586  );
1588  (
1589  objects,
1590  selectedFields
1591  );
1592  fvReconstructor.reconstructFvVolumeFields<symmTensor>
1593  (
1594  objects,
1595  selectedFields
1596  );
1597  fvReconstructor.reconstructFvVolumeFields<tensor>
1598  (
1599  objects,
1600  selectedFields
1601  );
1602 
1603 
1604  // surfaceFields
1605 
1606  fvReconstructor.reconstructFvSurfaceFields<scalar>
1607  (
1608  objects,
1609  selectedFields
1610  );
1611  fvReconstructor.reconstructFvSurfaceFields<vector>
1612  (
1613  objects,
1614  selectedFields
1615  );
1617  (
1618  objects,
1619  selectedFields
1620  );
1621  fvReconstructor.reconstructFvSurfaceFields<symmTensor>
1622  (
1623  objects,
1624  selectedFields
1625  );
1626  fvReconstructor.reconstructFvSurfaceFields<tensor>
1627  (
1628  objects,
1629  selectedFields
1630  );
1631 }
1633 
1635 (
1636  autoPtr<parLagrangianRedistributor>& lagrangianReconstructorPtr,
1637  const fvMesh& baseMesh,
1638  const fvMesh& mesh,
1639  const mapDistributePolyMesh& distMap,
1640  const HashSet<word>& selectedLagrangianFields
1641 )
1642 {
1643  // Clouds (note: might not be present on all processors)
1644 
1645  wordList cloudNames;
1647  // Find all cloudNames on all processors
1649 
1650  if (cloudNames.size())
1651  {
1652  if (!lagrangianReconstructorPtr.valid())
1653  {
1654  lagrangianReconstructorPtr.reset
1655  (
1657  (
1658  mesh,
1659  baseMesh,
1660  mesh.nCells(), // range of cell indices in clouds
1661  distMap
1662  )
1663  );
1664  }
1665  const parLagrangianRedistributor& lagrangianReconstructor =
1666  lagrangianReconstructorPtr();
1667 
1668  forAll(cloudNames, i)
1669  {
1670  Info<< "Reconstructing lagrangian fields for cloud "
1671  << cloudNames[i] << nl << endl;
1672 
1673  autoPtr<mapDistributeBase> lagrangianMap =
1674  lagrangianReconstructor.redistributeLagrangianPositions
1675  (
1676  cloudNames[i]
1677  );
1678  IOobjectList sprayObjs
1679  (
1680  mesh,
1681  mesh.time().timeName(),
1682  cloud::prefix/cloudNames[i]
1683  );
1684 
1685  lagrangianReconstructor.redistributeLagrangianFields<label>
1686  (
1687  lagrangianMap,
1688  cloudNames[i],
1689  sprayObjs,
1690  selectedLagrangianFields
1691  );
1692  lagrangianReconstructor.redistributeLagrangianFieldFields<label>
1693  (
1694  lagrangianMap,
1695  cloudNames[i],
1696  sprayObjs,
1697  selectedLagrangianFields
1698  );
1699  lagrangianReconstructor.redistributeLagrangianFields<scalar>
1700  (
1701  lagrangianMap,
1702  cloudNames[i],
1703  sprayObjs,
1704  selectedLagrangianFields
1705  );
1706  lagrangianReconstructor.redistributeLagrangianFieldFields<scalar>
1707  (
1708  lagrangianMap,
1709  cloudNames[i],
1710  sprayObjs,
1711  selectedLagrangianFields
1712  );
1713  lagrangianReconstructor.redistributeLagrangianFields<vector>
1714  (
1715  lagrangianMap,
1716  cloudNames[i],
1717  sprayObjs,
1718  selectedLagrangianFields
1719  );
1720  lagrangianReconstructor.redistributeLagrangianFieldFields<vector>
1721  (
1722  lagrangianMap,
1723  cloudNames[i],
1724  sprayObjs,
1725  selectedLagrangianFields
1726  );
1727  lagrangianReconstructor.redistributeLagrangianFields
1728  <sphericalTensor>
1729  (
1730  lagrangianMap,
1731  cloudNames[i],
1732  sprayObjs,
1733  selectedLagrangianFields
1734  );
1735  lagrangianReconstructor.redistributeLagrangianFieldFields
1736  <sphericalTensor>
1737  (
1738  lagrangianMap,
1739  cloudNames[i],
1740  sprayObjs,
1741  selectedLagrangianFields
1742  );
1743  lagrangianReconstructor.redistributeLagrangianFields<symmTensor>
1744  (
1745  lagrangianMap,
1746  cloudNames[i],
1747  sprayObjs,
1748  selectedLagrangianFields
1749  );
1750  lagrangianReconstructor.redistributeLagrangianFieldFields
1751  <symmTensor>
1752  (
1753  lagrangianMap,
1754  cloudNames[i],
1755  sprayObjs,
1756  selectedLagrangianFields
1757  );
1758  lagrangianReconstructor.redistributeLagrangianFields<tensor>
1759  (
1760  lagrangianMap,
1761  cloudNames[i],
1762  sprayObjs,
1763  selectedLagrangianFields
1764  );
1765  lagrangianReconstructor.redistributeLagrangianFieldFields<tensor>
1766  (
1767  lagrangianMap,
1768  cloudNames[i],
1769  sprayObjs,
1770  selectedLagrangianFields
1771  );
1772  }
1773  }
1774 }
1776 // Read clouds (note: might not be present on all processors)
1777 void readLagrangian
1778 (
1779  const fvMesh& mesh,
1780  const wordList& cloudNames,
1781  const HashSet<word>& selectedLagrangianFields,
1783 )
1784 {
1785  (void)mesh.tetBasePtIs();
1786 
1787  forAll(cloudNames, i)
1788  {
1789  {
1790  // Note: disable master-only reading of uniform/cloudProperties
1791  regIOobject::fileCheckTypes oldCheckType =
1793 
1794  if (oldCheckType == regIOobject::timeStampMaster)
1795  {
1797  }
1798  else if (oldCheckType == regIOobject::inotifyMaster)
1799  {
1801  }
1802 
1803  clouds.set
1804  (
1805  i,
1806  new unmappedPassiveParticleCloud(mesh, cloudNames[i], false)
1807  );
1808 
1810  }
1811 
1812 
1813  //forAllConstIter
1814  //(
1815  // unmappedPassiveParticleCloud,
1816  // clouds[i],
1817  // iter
1818  //)
1819  //{
1820  // Pout<< "Particle position:" << iter().position()
1821  // << " cell:" << iter().cell()
1822  // << " with cc:" << mesh.cellCentres()[iter().cell()]
1823  // << endl;
1824  //}
1825 
1826 
1827  IOobjectList sprayObjs(clouds[i], clouds[i].time().timeName());
1828 
1829  //Pout<< "Found clould objects:" << sprayObjs.names() << endl;
1830 
1832  <IOField<label> >
1833  (
1834  clouds[i],
1835  sprayObjs,
1836  selectedLagrangianFields
1837  );
1840  (
1841  clouds[i],
1842  sprayObjs,
1843  selectedLagrangianFields
1844  );
1847  (
1848  clouds[i],
1849  sprayObjs,
1850  selectedLagrangianFields
1851  );
1852 
1853 
1855  <IOField<scalar> >
1856  (
1857  clouds[i],
1858  sprayObjs,
1859  selectedLagrangianFields
1860  );
1863  (
1864  clouds[i],
1865  sprayObjs,
1866  selectedLagrangianFields
1867  );
1869  <CompactIOField<Field<scalar>, scalar> >
1870  (
1871  clouds[i],
1872  sprayObjs,
1873  selectedLagrangianFields
1874  );
1875 
1876 
1878  <IOField<vector> >
1879  (
1880  clouds[i],
1881  sprayObjs,
1882  selectedLagrangianFields
1883  );
1886  (
1887  clouds[i],
1888  sprayObjs,
1889  selectedLagrangianFields
1890  );
1893  (
1894  clouds[i],
1895  sprayObjs,
1896  selectedLagrangianFields
1897  );
1898 
1899 
1902  (
1903  clouds[i],
1904  sprayObjs,
1905  selectedLagrangianFields
1906  );
1909  (
1910  clouds[i],
1911  sprayObjs,
1912  selectedLagrangianFields
1913  );
1916  (
1917  clouds[i],
1918  sprayObjs,
1919  selectedLagrangianFields
1920  );
1921 
1922 
1925  (
1926  clouds[i],
1927  sprayObjs,
1928  selectedLagrangianFields
1929  );
1932  (
1933  clouds[i],
1934  sprayObjs,
1935  selectedLagrangianFields
1936  );
1939  (
1940  clouds[i],
1941  sprayObjs,
1942  selectedLagrangianFields
1943  );
1944 
1945 
1947  <IOField<tensor> >
1948  (
1949  clouds[i],
1950  sprayObjs,
1951  selectedLagrangianFields
1952  );
1955  (
1956  clouds[i],
1957  sprayObjs,
1958  selectedLagrangianFields
1959  );
1962  (
1963  clouds[i],
1964  sprayObjs,
1965  selectedLagrangianFields
1966  );
1967  }
1968 }
1970 
1972 (
1973  autoPtr<parLagrangianRedistributor>& lagrangianReconstructorPtr,
1974  const fvMesh& mesh,
1975  const label nOldCells,
1976  const mapDistributePolyMesh& distMap,
1978 )
1979 {
1980  if (clouds.size())
1981  {
1982  if (!lagrangianReconstructorPtr.valid())
1983  {
1984  lagrangianReconstructorPtr.reset
1985  (
1987  (
1988  mesh,
1989  mesh,
1990  nOldCells, // range of cell indices in clouds
1991  distMap
1992  )
1993  );
1994  }
1995  const parLagrangianRedistributor& distributor =
1996  lagrangianReconstructorPtr();
1997 
1998  forAll(clouds, i)
1999  {
2000  autoPtr<mapDistributeBase> lagrangianMap =
2001  distributor.redistributeLagrangianPositions(clouds[i]);
2002 
2004  <IOField<label> >
2005  (
2006  lagrangianMap,
2007  clouds[i]
2008  );
2011  (
2012  lagrangianMap,
2013  clouds[i]
2014  );
2017  (
2018  lagrangianMap,
2019  clouds[i]
2020  );
2021 
2022 
2024  <IOField<scalar> >
2025  (
2026  lagrangianMap,
2027  clouds[i]
2028  );
2031  (
2032  lagrangianMap,
2033  clouds[i]
2034  );
2036  <CompactIOField<Field<scalar>, scalar> >
2037  (
2038  lagrangianMap,
2039  clouds[i]
2040  );
2041 
2042 
2044  <IOField<vector> >
2045  (
2046  lagrangianMap,
2047  clouds[i]
2048  );
2051  (
2052  lagrangianMap,
2053  clouds[i]
2054  );
2057  (
2058  lagrangianMap,
2059  clouds[i]
2060  );
2061 
2062 
2065  (
2066  lagrangianMap,
2067  clouds[i]
2068  );
2071  (
2072  lagrangianMap,
2073  clouds[i]
2074  );
2077  (
2078  lagrangianMap,
2079  clouds[i]
2080  );
2081 
2082 
2085  (
2086  lagrangianMap,
2087  clouds[i]
2088  );
2091  (
2092  lagrangianMap,
2093  clouds[i]
2094  );
2097  (
2098  lagrangianMap,
2099  clouds[i]
2100  );
2101 
2102 
2104  <IOField<tensor> >
2105  (
2106  lagrangianMap,
2107  clouds[i]
2108  );
2111  (
2112  lagrangianMap,
2113  clouds[i]
2114  );
2117  (
2118  lagrangianMap,
2119  clouds[i]
2120  );
2121  }
2122  }
2124 
2125 
2126 int main(int argc, char *argv[])
2127 {
2128  // enable -constant ... if someone really wants it
2129  // enable -zeroTime to prevent accidentally trashing the initial fields
2130  timeSelector::addOptions(true, true);
2131  #include "addRegionOption.H"
2132  #include "addOverwriteOption.H"
2133  argList::addBoolOption("decompose", "decompose case");
2134  argList::addBoolOption("reconstruct", "reconstruct case");
2136  (
2137  "mergeTol",
2138  "scalar",
2139  "specify the merge distance relative to the bounding box size "
2140  "(default 1e-6)"
2141  );
2143  (
2144  "cellDist",
2145  "write cell distribution as a labelList - for use with 'manual' "
2146  "decomposition method or as a volScalarField for post-processing."
2147  );
2149  (
2150  "newTimes",
2151  "only reconstruct new times (i.e. that do not exist already)"
2152  );
2153 
2154 
2155  // Handle arguments
2156  // ~~~~~~~~~~~~~~~~
2157  // (replacement for setRootCase that does not abort)
2158 
2159  Foam::argList args(argc, argv);
2160  bool decompose = args.optionFound("decompose");
2161  const bool reconstruct = args.optionFound("reconstruct");
2162  bool overwrite = args.optionFound("overwrite");
2163  bool writeCellDist = args.optionFound("cellDist");
2164  bool newTimes = args.optionFound("newTimes");
2165 
2166 
2167 
2168  if (env("FOAM_SIGFPE"))
2169  {
2171  << "Detected floating point exception trapping (FOAM_SIGFPE)."
2172  << " This might give" << nl
2173  << " problems when mapping fields. Switch it off in case"
2174  << " of problems." << endl;
2175  }
2176 
2177 
2178 
2179  const HashSet<word> selectedFields(0);
2180  const HashSet<word> selectedLagrangianFields(0);
2181 
2182 
2183  if (decompose)
2184  {
2185  Info<< "Decomposing case (like decomposePar)" << nl << endl;
2186  if (reconstruct)
2187  {
2189  << "Cannot specify both -decompose and -reconstruct"
2190  << exit(FatalError);
2191  }
2192  }
2193  else if (reconstruct)
2194  {
2195  Info<< "Reconstructing case (like reconstructParMesh)" << nl << endl;
2196  }
2197 
2198 
2199  if (decompose || reconstruct)
2200  {
2201  if (!overwrite)
2202  {
2204  << "Working in decompose or reconstruction mode automatically"
2205  << " implies -overwrite" << nl << endl;
2206  overwrite = true;
2207  }
2208  }
2209 
2210 
2211  if (!Pstream::parRun())
2212  {
2214  << ": This utility can only be run parallel"
2215  << exit(FatalError);
2216  }
2217 
2218 
2219  if (!isDir(args.rootPath()))
2220  {
2222  << ": cannot open root directory " << args.rootPath()
2223  << exit(FatalError);
2224  }
2225 
2226  // Detect if running data-distributed (multiple roots)
2227  bool nfs = true;
2228  {
2229  List<fileName> roots(1, args.rootPath());
2231  nfs = (roots.size() == 1);
2232  }
2233 
2234  if (!nfs)
2235  {
2236  Info<< "Detected multiple roots i.e. non-nfs running"
2237  << nl << endl;
2238  }
2239 
2240  if (isDir(args.path()))
2241  {
2242  if (decompose)
2243  {
2244  Info<< "Removing existing processor directories" << endl;
2245  rmDir(args.path());
2246  }
2247  }
2248  else
2249  {
2250  // Directory does not exist. If this happens on master -> decompose mode
2251  decompose = true;
2252  Info<< "No processor directories; switching on decompose mode"
2253  << nl << endl;
2254  }
2255  // If master changed to decompose mode make sure all nodes know about
2256  // it
2257  Pstream::scatter(decompose);
2258 
2259 
2260 
2261 
2262  // If running distributed we have problem of new processors not finding
2263  // a system/controlDict. However if we switch on the master-only reading
2264  // the problem becomes that the time directories are differing sizes and
2265  // e.g. latestTime will pick up a different time (which causes createTime.H
2266  // to abort). So for now make sure to have master times on all
2267  // processors
2268  {
2269  Info<< "Creating time directories on all processors" << nl << endl;
2271  if (Pstream::master())
2272  {
2273  timeDirs = Time::findTimes(args.path(), "constant");
2274  }
2276  forAll(timeDirs, i)
2277  {
2278  mkDir(args.path()/timeDirs[i].name());
2279  }
2280  }
2281 
2282 
2283  // Construct time
2284  // ~~~~~~~~~~~~~~
2285 
2286  #include "createTime.H"
2287  runTime.functionObjects().off();
2288 
2289 
2290  // Save local processor0 casename
2291  const fileName proc0CaseName = runTime.caseName();
2292 
2293 
2294  // Construct undecomposed Time
2295  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
2296  // This will read the same controlDict but might have a different
2297  // set of times so enforce same times
2298 
2299  if (!nfs)
2300  {
2301  Info<< "Creating time directories for undecomposed Time"
2302  << " on all processors" << nl << endl;
2304 
2305  const fileName basePath(args.rootPath()/args.globalCaseName());
2306 
2307  if (Pstream::master())
2308  {
2309  timeDirs = Time::findTimes(basePath, "constant");
2310  }
2312  forAll(timeDirs, i)
2313  {
2314  mkDir(basePath/timeDirs[i].name());
2315  }
2316  }
2317 
2318 
2319  Info<< "Create undecomposed database"<< nl << endl;
2320  Time baseRunTime
2321  (
2322  runTime.controlDict(),
2323  runTime.rootPath(),
2324  runTime.globalCaseName(),
2325  runTime.system(),
2326  runTime.constant(),
2327  false // enableFunctionObjects
2328  );
2329 
2330 
2331  HashSet<word> masterTimeDirSet;
2332  if (newTimes)
2333  {
2334  instantList baseTimeDirs(baseRunTime.times());
2335  forAll(baseTimeDirs, i)
2336  {
2337  masterTimeDirSet.insert(baseTimeDirs[i].name());
2338  }
2339  }
2340 
2341 
2342  // Determine any region
2344  fileName meshSubDir;
2345 
2346  if (args.optionReadIfPresent("region", regionName))
2347  {
2348  meshSubDir = regionName/polyMesh::meshSubDir;
2349  }
2350  else
2351  {
2352  meshSubDir = polyMesh::meshSubDir;
2353  }
2354  Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
2355 
2356 
2357 
2358  // Demand driven lagrangian mapper
2359  autoPtr<parLagrangianRedistributor> lagrangianReconstructorPtr;
2360 
2361 
2362 
2363  if (reconstruct)
2364  {
2365  // use the times list from the master processor
2366  // and select a subset based on the command-line options
2369 
2370  if (timeDirs.empty())
2371  {
2373  << "No times selected"
2374  << exit(FatalError);
2375  }
2376 
2377 
2378  // Pass1 : reconstruct mesh and addressing
2379  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2380 
2381 
2382  Info<< nl
2383  << "Pass1 : reconstructing mesh and addressing" << nl << endl;
2384 
2385  // Loop over all times
2386  forAll(timeDirs, timeI)
2387  {
2388  // Set time for global database
2389  runTime.setTime(timeDirs[timeI], timeI);
2390  baseRunTime.setTime(timeDirs[timeI], timeI);
2391 
2392  Info<< "Time = " << runTime.timeName() << endl << endl;
2393 
2394 
2395  // See where the mesh is
2396  fileName facesInstance = runTime.findInstance
2397  (
2398  meshSubDir,
2399  "faces",
2401  );
2402  //Pout<< "facesInstance:" << facesInstance << endl;
2403 
2404  Pstream::scatter(facesInstance);
2405 
2406  // Check who has a mesh
2407  const fileName meshPath =
2408  runTime.path()/facesInstance/meshSubDir/"faces";
2409 
2410  Info<< "Checking for mesh in " << meshPath << nl << endl;
2411 
2412  boolList haveMesh(Pstream::nProcs(), false);
2413  haveMesh[Pstream::myProcNo()] = isFile(meshPath);
2414  Pstream::gatherList(haveMesh);
2415  Pstream::scatterList(haveMesh);
2416  Info<< "Per processor mesh availability : " << haveMesh << endl;
2417 
2418 
2419  // Addressing back to reconstructed mesh as xxxProcAddressing.
2420  // - all processors have consistent faceProcAddressing
2421  // - processors with no mesh don't need faceProcAddressing
2422 
2423 
2424  // Note: filePath searches up on processors that don't have
2425  // processor if instance = constant so explicitly check found
2426  // filename.
2427  bool haveAddressing = false;
2428  if (haveMesh[Pstream::myProcNo()])
2429  {
2430  haveAddressing = IOobject
2431  (
2432  "faceProcAddressing",
2433  facesInstance,
2434  meshSubDir,
2435  runTime,
2437  ).headerOk();
2438  }
2439  else
2440  {
2441  // Have no mesh. Don't need addressing
2442  haveAddressing = true;
2443  }
2444 
2445  if (!returnReduce(haveAddressing, andOp<bool>()))
2446  {
2447  Info<< "loading mesh from " << facesInstance << endl;
2449  (
2450  IOobject
2451  (
2452  regionName,
2453  facesInstance,
2454  runTime,
2456  )
2457  );
2458  fvMesh& mesh = meshPtr();
2459 
2460  // Global matching tolerance
2461  const scalar tolDim = getMergeDistance
2462  (
2463  args,
2464  runTime,
2465  mesh.bounds()
2466  );
2467 
2468 
2469  // Determine decomposition
2470  // ~~~~~~~~~~~~~~~~~~~~~~~
2471 
2472  Info<< "Reconstructing mesh for time " << facesInstance << endl;
2473 
2474  label nDestProcs = 1;
2475  labelList finalDecomp = labelList(mesh.nCells(), 0);
2476 
2478  (
2479  baseRunTime,
2480  tolDim,
2481  haveMesh,
2482  meshSubDir,
2483  false, // do not read fields
2484  false, // do not read undecomposed case on processor0
2485  overwrite,
2486  proc0CaseName,
2487  nDestProcs,
2488  finalDecomp,
2489  facesInstance,
2490  mesh
2491  );
2492  }
2493  }
2494 
2495 
2496  // Pass2 : read mesh and addressing and reconstruct fields
2497  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2498 
2499  Info<< nl
2500  << "Pass2 : reconstructing fields" << nl << endl;
2501 
2502  runTime.setTime(timeDirs[0], 0);
2503  baseRunTime.setTime(timeDirs[0], 0);
2504  Info<< "Time = " << runTime.timeName() << endl << endl;
2505 
2506 
2507  Info<< "Reading undecomposed mesh (on master)" << endl;
2509  (
2510  IOobject
2511  (
2512  regionName,
2513  baseRunTime.timeName(),
2514  baseRunTime,
2516  ),
2517  true // read on master only
2518  );
2519 
2520  Info<< "Reading local, decomposed mesh" << endl;
2522  (
2523  IOobject
2524  (
2525  regionName,
2526  baseMeshPtr().facesInstance(),
2527  runTime,
2529  )
2530  );
2531  fvMesh& mesh = meshPtr();
2532 
2533 
2534  // Read addressing back to base mesh
2536  readProcAddressing(mesh, baseMeshPtr, distMap);
2537 
2538  // Construct field mapper
2539  autoPtr<parFvFieldReconstructor> fvReconstructorPtr
2540  (
2542  (
2543  baseMeshPtr(),
2544  mesh,
2545  distMap,
2546  Pstream::master() // do I need to write?
2547  )
2548  );
2549 
2550 
2551 
2552  // Since we start from Times[0] and not runTime.timeName() we
2553  // might overlook point motion in the first timestep
2554  // (since mesh.readUpdate() below will not be triggered). Instead
2555  // detect points by hand
2557  {
2558  Info<< " Dected initial mesh motion; reconstructing points" << nl
2559  << endl;
2560  fvReconstructorPtr().reconstructPoints();
2561  }
2562 
2563 
2564  // Loop over all times
2565  forAll(timeDirs, timeI)
2566  {
2567  if (newTimes && masterTimeDirSet.found(timeDirs[timeI].name()))
2568  {
2569  Info<< "Skipping time " << timeDirs[timeI].name()
2570  << endl << endl;
2571  continue;
2572  }
2573 
2574  // Set time for global database
2575  runTime.setTime(timeDirs[timeI], timeI);
2576  baseRunTime.setTime(timeDirs[timeI], timeI);
2577 
2578  Info<< "Time = " << runTime.timeName() << endl << endl;
2579 
2580 
2581  // Check if any new meshes need to be read.
2583 
2584  if (procStat == fvMesh::POINTS_MOVED)
2585  {
2586  Info<< " Dected mesh motion; reconstructing points" << nl
2587  << endl;
2588  fvReconstructorPtr().reconstructPoints();
2589  }
2590  else if
2591  (
2592  procStat == fvMesh::TOPO_CHANGE
2593  || procStat == fvMesh::TOPO_PATCH_CHANGE
2594  )
2595  {
2596  Info<< " Detected topology change; reconstructing addressing"
2597  << nl << endl;
2598 
2599  if (baseMeshPtr.valid())
2600  {
2601  // Cannot do a baseMesh::readUpdate() since not all
2602  // processors will have mesh files. So instead just
2603  // recreate baseMesh
2604  baseMeshPtr.clear();
2605  baseMeshPtr = fvMeshTools::newMesh
2606  (
2607  IOobject
2608  (
2609  regionName,
2610  baseRunTime.timeName(),
2611  baseRunTime,
2613  ),
2614  true // read on master only
2615  );
2616  }
2617 
2618  // Re-read procXXXaddressing
2619  readProcAddressing(mesh, baseMeshPtr, distMap);
2620 
2621  // Reset field mapper
2622  fvReconstructorPtr.reset
2623  (
2625  (
2626  baseMeshPtr(),
2627  mesh,
2628  distMap,
2629  Pstream::master()
2630  )
2631  );
2632  lagrangianReconstructorPtr.clear();
2633  }
2634 
2635 
2636  // Get list of objects
2637  IOobjectList objects(mesh, runTime.timeName());
2638 
2639 
2640  // Mesh fields (vol, surface, volInternal)
2642  (
2643  fvReconstructorPtr(),
2644  objects,
2645  selectedFields
2646  );
2647 
2648  // Clouds (note: might not be present on all processors)
2650  (
2651  lagrangianReconstructorPtr,
2652  baseMeshPtr(),
2653  mesh,
2654  distMap,
2655  selectedLagrangianFields
2656  );
2657 
2658  // If there are any "uniform" directories copy them from
2659  // the master processor
2660  if (Pstream::master())
2661  {
2662  fileName uniformDir0 = runTime.timePath()/"uniform";
2663  if (isDir(uniformDir0))
2664  {
2665  Info<< "Detected additional non-decomposed files in "
2666  << uniformDir0 << endl;
2667  cp(uniformDir0, baseRunTime.timePath());
2668  }
2669  }
2670  }
2671  }
2672  else
2673  {
2674  // Time coming from processor0 (or undecomposed if no processor0)
2675  scalar masterTime;
2676  if (decompose)
2677  {
2678  // Use base time. This is to handle e.g. startTime = latestTime
2679  // which will not do anything if there are no processor directories
2680  masterTime = timeSelector::selectIfPresent
2681  (
2682  baseRunTime,
2683  args
2684  )[0].value();
2685  }
2686  else
2687  {
2688  masterTime = timeSelector::selectIfPresent
2689  (
2690  runTime,
2691  args
2692  )[0].value();
2693  }
2694  Pstream::scatter(masterTime);
2695  Info<< "Setting time to that of master or undecomposed case : "
2696  << masterTime << endl;
2697  runTime.setTime(masterTime, 0);
2698 
2699 
2700  // Get time instance directory
2701  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
2702  // At this point we should be able to read at least a mesh on
2703  // processor0. Note the changing of the processor0 casename to
2704  // enforce it to read/write from the undecomposed case
2705 
2706  fileName masterInstDir;
2707  if (Pstream::master())
2708  {
2709  if (decompose)
2710  {
2711  Info<< "Setting caseName to " << baseRunTime.caseName()
2712  << " to find undecomposed mesh" << endl;
2713  runTime.TimePaths::caseName() = baseRunTime.caseName();
2714  }
2715 
2716  masterInstDir = runTime.findInstance
2717  (
2718  meshSubDir,
2719  "faces",
2721  );
2722 
2723  if (decompose)
2724  {
2725  Info<< "Restoring caseName to " << proc0CaseName << endl;
2726  runTime.TimePaths::caseName() = proc0CaseName;
2727  }
2728  }
2729  Pstream::scatter(masterInstDir);
2730 
2731  // Check who has a mesh
2732  const fileName meshPath =
2733  runTime.path()/masterInstDir/meshSubDir/"faces";
2734 
2735  Info<< "Checking for mesh in " << meshPath << nl << endl;
2736 
2737 
2738  boolList haveMesh(Pstream::nProcs(), false);
2739  haveMesh[Pstream::myProcNo()] = isFile(meshPath);
2740  Pstream::gatherList(haveMesh);
2741  Pstream::scatterList(haveMesh);
2742  Info<< "Per processor mesh availability : " << haveMesh << endl;
2743 
2744  // Load mesh (or create dummy one)
2745  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2746 
2747  if (Pstream::master() && decompose)
2748  {
2749  Info<< "Setting caseName to " << baseRunTime.caseName()
2750  << " to read undecomposed mesh" << endl;
2751  runTime.TimePaths::caseName() = baseRunTime.caseName();
2752  }
2753 
2755  (
2756  IOobject
2757  (
2758  regionName,
2759  masterInstDir,
2760  runTime,
2762  )
2763  );
2764 
2765  if (Pstream::master() && decompose)
2766  {
2767  Info<< "Restoring caseName to " << proc0CaseName << endl;
2768  runTime.TimePaths::caseName() = proc0CaseName;
2769  }
2770 
2771  fvMesh& mesh = meshPtr();
2772 
2773 
2774  label nOldCells = mesh.nCells();
2775  //Pout<< "Loaded mesh : nCells:" << nOldCells
2776  // << " nPatches:" << mesh.boundaryMesh().size() << endl;
2777 
2778 
2779  // Global matching tolerance
2780  const scalar tolDim = getMergeDistance
2781  (
2782  args,
2783  runTime,
2784  mesh.bounds()
2785  );
2786 
2787  // Allow override of decomposeParDict location
2788  fileName decompDictFile;
2789  if (args.optionReadIfPresent("decomposeParDict", decompDictFile))
2790  {
2791  if (isDir(decompDictFile))
2792  {
2793  decompDictFile = decompDictFile / "decomposeParDict";
2794  }
2795  }
2796 
2797 
2798  // Determine decomposition
2799  // ~~~~~~~~~~~~~~~~~~~~~~~
2800 
2801  label nDestProcs;
2802  labelList finalDecomp;
2804  (
2805  baseRunTime,
2806  decompDictFile,
2807  decompose,
2808  proc0CaseName,
2809  mesh,
2810  writeCellDist,
2811 
2812  nDestProcs,
2813  finalDecomp
2814  );
2815 
2816 
2817  wordList cloudNames;
2819 
2820  // Detect lagrangian fields
2821  if (Pstream::master() && decompose)
2822  {
2823  runTime.TimePaths::caseName() = baseRunTime.caseName();
2824  }
2826  (
2827  mesh,
2828  cloudNames,
2829  fieldNames
2830  );
2831 
2832  // Read lagrangian fields and store on cloud (objectRegistry)
2833  PtrList<unmappedPassiveParticleCloud> clouds(cloudNames.size());
2835  (
2836  mesh,
2837  cloudNames,
2838  selectedLagrangianFields,
2839  clouds
2840  );
2841  if (Pstream::master() && decompose)
2842  {
2843  runTime.TimePaths::caseName() = proc0CaseName;
2844  }
2845 
2846 
2847  // Load fields, do all distribution (mesh and fields - but not
2848  // lagrangian fields; these are done later)
2850  (
2851  baseRunTime,
2852  tolDim,
2853  haveMesh,
2854  meshSubDir,
2855  true, // read fields
2856  decompose, // decompose, i.e. read from undecomposed case
2857  overwrite,
2858  proc0CaseName,
2859  nDestProcs,
2860  finalDecomp,
2861  masterInstDir,
2862  mesh
2863  );
2864 
2865 
2866  // Redistribute any clouds
2868  (
2869  lagrangianReconstructorPtr,
2870  mesh,
2871  nOldCells,
2872  distMap,
2873  clouds
2874  );
2875 
2876 
2877  // Copy any uniform data
2878  const fileName uniformDir("uniform");
2879  if (isDir(baseRunTime.timePath()/uniformDir))
2880  {
2881  Info<< "Detected additional non-decomposed files in "
2882  << baseRunTime.timePath()/uniformDir << endl;
2883  cp
2884  (
2885  baseRunTime.timePath()/uniformDir,
2886  runTime.timePath()/uniformDir
2887  );
2888  }
2889  }
2890 
2891 
2892  Info<< "End\n" << endl;
2893 
2894  return 0;
2895 }
2896 
2897 
2898 // ************************************************************************* //
createReconstructMap
autoPtr< mapDistributePolyMesh > createReconstructMap(const autoPtr< fvMesh > &baseMeshPtr, const fvMesh &mesh, const labelList &cellProcAddressing, const labelList &faceProcAddressing, const labelList &pointProcAddressing, const labelList &boundaryProcAddressing)
Definition: redistributePar.C:1220
volFields.H
Foam::unmappedPassiveParticleCloud
passiveParticleCloud but with autoMap and writing disabled. Only used for its objectRegistry to store...
Definition: unmappedPassiveParticleCloud.H:49
Foam::mapDistributeBase::subMap
const labelListList & subMap() const
From subsetted data back to original data.
Definition: mapDistributeBase.H:256
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:54
reconstructMeshFields
void reconstructMeshFields(const parFvFieldReconstructor &fvReconstructor, const IOobjectList &objects, const HashSet< word > &selectedFields)
Definition: redistributePar.C:1537
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::mapDistributeBase::constructHasFlip
bool constructHasFlip() const
Does constructMap include a sign.
Definition: mapDistributeBase.H:292
Foam::parLagrangianRedistributor::redistributeStoredLagrangianFields
void redistributeStoredLagrangianFields(const mapDistributeBase &map, passiveParticleCloud &cloud) const
Redistribute and write stored lagrangian fields.
Definition: parLagrangianRedistributorRedistributeFields.C:273
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::env
bool env(const word &)
Return true if environment variable of given name is defined.
Definition: POSIX.C:95
Foam::SymmTensor< scalar >
Foam::boundBox::mag
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
Foam::parLagrangianRedistributor::readLagrangianFields
static void readLagrangianFields(const passiveParticleCloud &cloud, const IOobjectList &objects, const HashSet< word > &selectedFields)
Read and store all fields of a cloud.
Definition: parLagrangianRedistributorRedistributeFields.C:225
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
Foam::parLagrangianRedistributor::redistributeLagrangianPositions
autoPtr< mapDistributeBase > redistributeLagrangianPositions(passiveParticleCloud &cloud) const
Redistribute and write lagrangian positions.
Definition: parLagrangianRedistributor.C:130
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
processorFvPatchField.H
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
loadOrCreateMesh.H
Load or create (0 size) a mesh. Used in distributing meshes to a larger number of processors.
Foam::polyMesh::POINTS_MOVED
@ POINTS_MOVED
Definition: polyMesh.H:91
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::argList::addOption
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:108
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Foam::UPstream::scheduled
@ scheduled
Definition: UPstream.H:67
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::processorFvPatchField
This boundary condition enables processor communication across patches.
Definition: processorFvPatchField.H:64
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
Foam::Time::caseName
const fileName & caseName() const
Return case name.
Definition: Time.H:275
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
readField
void readField(const IOobject &io, const fvMesh &mesh, const label i, PtrList< GeoField > &fields)
Definition: redistributePar.C:557
Foam::DynamicList< label >
Foam::globalMeshData::processorPatches
const labelList & processorPatches() const
Return list of processor patch labels.
Definition: globalMeshData.H:404
Foam::mapDistributeBase::distribute
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &, const negateOp &negOp, const int tag=UPstream::msgType())
Distribute data. Note:schedule only used for Pstream::scheduled.
Definition: mapDistributeBaseTemplates.C:119
Foam::List::xfer
Xfer< List< T > > xfer()
Transfer contents to the Xfer container.
Definition: ListI.H:90
writeDecomposition
void writeDecomposition(const word &name, const fvMesh &mesh, const labelList &decomp)
Definition: redistributePar.C:226
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:50
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
Foam::Time::times
instantList times() const
Search the case for valid time directories.
Definition: Time.C:758
Foam::polyMesh::tetBasePtIs
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:818
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:205
Foam::Time::writeFormat
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:303
Foam::mapDistributePolyMesh::distributeCellData
void distributeCellData(List< T > &lst) const
Distribute list of cell data.
Definition: mapDistributePolyMesh.H:250
addOverwriteOption.H
correctCoupledBoundaryConditions
void correctCoupledBoundaryConditions(fvMesh &mesh)
Definition: redistributePar.C:700
globalIndex.H
IOmapDistributePolyMesh.H
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::regIOobject::inotify
@ inotify
Definition: regIOobject.H:72
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::flipOp
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:50
printMeshData
void printMeshData(const polyMesh &mesh)
Definition: redistributePar.C:128
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
Foam::UPstream::commsTypeNames
static const NamedEnum< commsTypes, 3 > commsTypeNames
Definition: UPstream.H:71
Foam::Time::functionObjects
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:435
Foam::IOobjectList::sortedNames
wordList sortedNames() const
Return the sorted list of names of the IOobjects.
Definition: IOobjectList.C:227
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::UPstream::waitRequests
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:106
Foam::combineReduce
void combineReduce(const List< UPstream::commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
Definition: PstreamCombineReduceOps.H:52
Foam::processorPolyPatch::neighbProcNo
int neighbProcNo() const
Return neigbour processor number.
Definition: processorPolyPatch.H:255
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:147
Foam::argList
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:97
Foam::globalIndex::localSize
label localSize() const
My local size.
Definition: globalIndexI.H:60
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
IOobjectList.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::UPstream::defaultCommsType
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:252
Foam::CompactIOField
A Field of objects of type <T> with automated input and output using a compact storage....
Definition: CompactIOField.H:50
decompositionMethod.H
Foam::IOobjectList::names
wordList names() const
Return the list of names of the IOobjects.
Definition: IOobjectList.C:221
Foam::cp
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:755
Foam::IOobject::writeOpt
writeOption writeOpt() const
Definition: IOobject.H:327
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
Foam::decompositionModel
MeshObject wrapper of decompositionMethod.
Definition: decompositionModel.H:53
Foam::argList::rootPath
const fileName & rootPath() const
Return root path.
Definition: argListI.H:36
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
parFvFieldReconstructor.H
Foam::mapDistributePolyMesh::cellMap
const mapDistribute & cellMap() const
Cell distribute map.
Definition: mapDistributePolyMesh.H:214
Foam::parLagrangianRedistributor::redistributeLagrangianFields
void redistributeLagrangianFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, redistribute and write all/selected lagrangian fields.
Definition: parLagrangianRedistributorRedistributeFields.C:68
Foam::fvMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:498
Foam::IOobject::objectPath
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:376
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
getMergeDistance
scalar getMergeDistance(const argList &args, const Time &runTime, const boundBox &bb)
Definition: redistributePar.C:89
Foam::polyMesh::TOPO_PATCH_CHANGE
@ TOPO_PATCH_CHANGE
Definition: polyMesh.H:93
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::regIOobject::timeStamp
@ timeStamp
Definition: regIOobject.H:70
Foam::loadOrCreateMesh
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
Definition: loadOrCreateMesh.C:46
Foam::decompositionMethod::parallelAware
virtual bool parallelAware() const =0
Is method parallel aware (i.e. does it synchronize domains across.
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::IOstream::ASCII
@ ASCII
Definition: IOstream.H:88
Foam::Time::rootPath
const fileName & rootPath() const
Return root path.
Definition: Time.H:269
Foam::PtrList::set
bool set(const label) const
Is element set.
determineDecomposition
void determineDecomposition(const Time &baseRunTime, const fileName &decompDictFile, const bool decompose, const fileName &proc0CaseName, const fvMesh &mesh, const bool writeCellDist, label &nDestProcs, labelList &decomp)
Definition: redistributePar.C:277
redistributeLagrangian
void redistributeLagrangian(autoPtr< parLagrangianRedistributor > &lagrangianReconstructorPtr, const fvMesh &mesh, const label nOldCells, const mapDistributePolyMesh &distMap, PtrList< unmappedPassiveParticleCloud > &clouds)
Definition: redistributePar.C:1969
Foam::UPstream::nonBlocking
@ nonBlocking
Definition: UPstream.H:68
Foam::label
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
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::TimePaths::system
const word & system() const
Return system name.
Definition: TimePaths.H:120
Foam::functionObjectList::off
virtual void off()
Switch the function objects off.
Definition: functionObjectList.C:178
Foam::parFvFieldReconstructor::reconstructFvVolumeInternalFields
void reconstructFvVolumeInternalFields(const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, reconstruct and write all/selected volume internal fields.
Definition: parFvFieldReconstructorReconstructFields.C:420
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::parLagrangianRedistributor
Lagrangian field redistributor.
Definition: parLagrangianRedistributor.H:57
Foam::Time::timePath
fileName timePath() const
Return current time path.
Definition: Time.H:297
Foam::fvMeshTools::newMesh
static autoPtr< fvMesh > newMesh(const IOobject &io, const bool masterOnlyReading)
Read mesh or create dummy mesh (0 cells, >0 patches). Works in two.
Definition: fvMeshTools.C:422
parLagrangianRedistributor.H
Foam::mapDistributePolyMesh::faceMap
const mapDistribute & faceMap() const
Face distribute map.
Definition: mapDistributePolyMesh.H:208
argList.H
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
Foam::decompositionModel::decomposer
decompositionMethod & decomposer() const
Definition: decompositionModel.H:110
Foam::mapDistributePolyMesh::patchMap
const mapDistribute & patchMap() const
Patch distribute map.
Definition: mapDistributePolyMesh.H:220
addRegionOption.H
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:55
Foam::flipLabelOp
Definition: flipOp.H:74
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::polyMesh::TOPO_CHANGE
@ TOPO_CHANGE
Definition: polyMesh.H:92
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::UPstream::masterNo
static int masterNo()
Process index of the master.
Definition: UPstream.H:393
Foam::hexRef8Data::distribute
void distribute(const mapDistributePolyMesh &)
In-place distribute.
Definition: hexRef4Data.C:304
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
readLagrangian
void readLagrangian(const fvMesh &mesh, const wordList &cloudNames, const HashSet< word > &selectedLagrangianFields, PtrList< unmappedPassiveParticleCloud > &clouds)
Definition: redistributePar.C:1775
Foam::hexRef8Data
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition: hexRef4Data.H:57
Foam::regIOobject::fileCheckTypes
fileCheckTypes
Types of communications.
Definition: regIOobject.H:68
Foam::timeSelector::selectIfPresent
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times (as select0)
Definition: timeSelector.C:284
Foam::mapDistributePolyMesh::oldPatchSizes
const labelList & oldPatchSizes() const
List of the old patch sizes.
Definition: mapDistributePolyMesh.H:184
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Foam::argList::path
fileName path() const
Return the path to the caseName.
Definition: argListI.H:60
meshPtr
static fvMesh * meshPtr
Definition: globalFoam.H:52
Foam::parLagrangianRedistributor::findClouds
static void findClouds(const fvMesh &, wordList &cloudNames, List< wordList > &objectNames)
Find all clouds (on all processors) and for each cloud all.
Definition: parLagrangianRedistributor.C:60
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::decompositionMethod::nDomains
label nDomains() const
Definition: decompositionMethod.H:112
Foam::Time::controlDict
const dictionary & controlDict() const
Definition: Time.H:286
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
readProcAddressing
void readProcAddressing(const fvMesh &mesh, const autoPtr< fvMesh > &baseMeshPtr, autoPtr< mapDistributePolyMesh > &distMap)
Definition: redistributePar.C:1422
Foam::Time::setTime
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:964
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::mapDistributePolyMesh::xfer
Xfer< mapDistributePolyMesh > xfer()
Transfer contents to the Xfer container.
Definition: mapDistributePolyMesh.C:197
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::decompositionMethod
Abstract base class for decomposition.
Definition: decompositionMethod.H:48
Foam::eqOp
Definition: ops.H:70
Foam::mapDistributeBase::subHasFlip
bool subHasFlip() const
Does subMap include a sign.
Definition: mapDistributeBase.H:280
defaultMergeTol
static const scalar defaultMergeTol
Definition: redistributePar.C:84
Foam::HashTable< nil, word, string::hash >::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
hexRef8Data.H
Foam::mapDistributePolyMesh::transfer
void transfer(mapDistributePolyMesh &)
Transfer the contents of the argument and annul the argument.
Definition: mapDistributePolyMesh.C:182
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::isFile
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:622
Foam::isDir
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:615
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:199
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::parFvFieldReconstructor::reconstructFvSurfaceFields
void reconstructFvSurfaceFields(const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, reconstruct and write all/selected surface fields.
Definition: parFvFieldReconstructorReconstructFields.C:504
PstreamReduceOps.H
Foam::SphericalTensor
Templated 3D SphericalTensor derived from VectorSpace adding construction from 1 component,...
Definition: SphericalTensor.H:51
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Foam::Time::findInstance
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Return the location of "dir" containing the file "name".
Definition: findInstance.C:38
Foam::hexRef8Data::sync
void sync(const IOobject &io)
Parallel synchronise. This enforces valid objects on all processors.
Definition: hexRef4Data.C:246
Foam::ListUniqueEqOp
Helper class for list to append unique elelements of y onto the end of x.
Definition: ListOps.H:268
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
Foam::cloud::prefix
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
Foam::regIOobject::fileModificationChecking
static fileCheckTypes fileModificationChecking
Definition: regIOobject.H:128
Foam::mapDistributeBase::constructSize
label constructSize() const
Constructed data size.
Definition: mapDistributeBase.H:244
Foam::HashPtrTable< IOobject >::iterator
HashTable< T *, Key, Hash >::iterator iterator
Definition: HashPtrTable.H:82
Foam::decompositionModel::New
static const decompositionModel & New(const polyMesh &mesh, const fileName &decompDictFile="")
Read (optionallly from absolute path) & register on mesh.
Definition: decompositionModel.C:108
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::andOp
Definition: ops.H:177
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
Foam::decompositionMethod::decompose
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Return for every coordinate the wanted processor number.
Definition: decompositionMethod.H:126
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:452
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:49
Foam::hexRef8Data::write
bool write() const
Write.
Definition: hexRef4Data.C:324
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::UPstream::nRequests
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:96
fvMeshDistribute.H
Foam::regIOobject::inotifyMaster
@ inotifyMaster
Definition: regIOobject.H:73
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:281
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
fvMeshTools.H
Foam::TimePaths::caseName
const fileName & caseName() const
Return case name.
Definition: TimePaths.H:108
Foam::globalIndex::size
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Foam::Vector< scalar >
Foam::parLagrangianRedistributor::redistributeLagrangianFieldFields
void redistributeLagrangianFieldFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, redistribute and write all/selected lagrangian fieldFields.
Definition: parLagrangianRedistributorRedistributeFields.C:140
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::fieldNames
wordList fieldNames(const IOobjectList &objects, const bool syncPar)
Get sorted names of fields of type. If syncPar and running in parallel.
Definition: ReadFields.C:36
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::IOmapDistributePolyMesh
IOmapDistributePolyMesh is derived from mapDistributePolyMesh and IOobject to give the mapDistributeP...
Definition: IOmapDistributePolyMesh.H:51
Foam::Time::findTimes
static instantList findTimes(const fileName &, const word &constantName="constant")
Search a given directory for valid time directories.
Definition: findTimes.C:38
main
int main(int argc, char *argv[])
Definition: redistributePar.C:2123
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePaths.H:130
Foam::autoPtr::valid
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set).
Definition: autoPtrI.H:83
Foam::HashPtrTable::erase
bool erase(iterator &)
Erase an hashedEntry specified by given iterator.
Definition: HashPtrTable.C:76
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::IOList< label >
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
timeName
word timeName
Definition: getTimeIndex.H:3
timeSelector.H
createTime.H
Foam::fvMeshDistribute::distribute
autoPtr< mapDistributePolyMesh > distribute(const labelList &dist)
Send cells to neighbours according to distribution.
Definition: fvMeshDistribute.C:1613
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
faceProcAddressing
PtrList< labelIOList > & faceProcAddressing
Definition: checkFaceAddressingComp.H:9
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::autoPtr::clear
void clear()
Delete object (if the pointer is valid) and set pointer to NULL.
Definition: autoPtrI.H:126
writeProcAddressing
void writeProcAddressing(const bool decompose, const fileName &meshSubDir, const fvMesh &mesh, const mapDistributePolyMesh &map)
Definition: redistributePar.C:378
Foam::mapDistribute::reverseDistribute
void reverseDistribute(const label constructSize, List< T > &, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
Definition: mapDistributeTemplates.C:187
Foam::autoPtr::reset
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Foam::timeSelector::select
instantList select(const instantList &) const
Select a list of Time values that are within the ranges.
Definition: timeSelector.C:100
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
unmappedPassiveParticleCloud.H
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
Foam::mapDistributePolyMesh::pointMap
const mapDistribute & pointMap() const
Point distribute map.
Definition: mapDistributePolyMesh.H:202
decompositionModel.H
reconstructLagrangian
void reconstructLagrangian(autoPtr< parLagrangianRedistributor > &lagrangianReconstructorPtr, const fvMesh &baseMesh, const fvMesh &mesh, const mapDistributePolyMesh &distMap, const HashSet< word > &selectedLagrangianFields)
Definition: redistributePar.C:1632
Foam::mapDistributePolyMesh::nOldFaces
label nOldFaces() const
Number of faces in mesh before distribution.
Definition: mapDistributePolyMesh.H:172
Foam::volSphericalTensorField
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
Definition: volFieldsFwd.H:57
Foam::parFvFieldReconstructor::reconstructFvVolumeFields
void reconstructFvVolumeFields(const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, reconstruct and write all/selected volume fields.
Definition: parFvFieldReconstructorReconstructFields.C:461
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::TimePaths::globalCaseName
const fileName & globalCaseName() const
Return global case name.
Definition: TimePaths.H:102
Foam::mapDistributePolyMesh
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Definition: mapDistributePolyMesh.H:57
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::mapDistributePolyMesh::nOldPoints
label nOldPoints() const
Number of points in mesh before distribution.
Definition: mapDistributePolyMesh.H:166
Foam::mapDistributePolyMesh::nOldCells
label nOldCells() const
Number of cells in mesh before distribution.
Definition: mapDistributePolyMesh.H:178
redistributeAndWrite
autoPtr< mapDistributePolyMesh > redistributeAndWrite(const Time &baseRunTime, const scalar tolDim, const boolList &haveMesh, const fileName &meshSubDir, const bool doReadFields, const bool decompose, const bool overwrite, const fileName &proc0CaseName, const label nDestProcs, const labelList &decomp, const fileName &masterInstDir, fvMesh &mesh)
Definition: redistributePar.C:793
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::mapDistributeBase::constructMap
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
Definition: mapDistributeBase.H:268
Foam::mapDistribute::xfer
Xfer< mapDistribute > xfer()
Transfer contents to the Xfer container.
Definition: mapDistribute.C:531
args
Foam::argList args(argc, argv)
Foam::mkDir
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:419
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::mapDistributePolyMesh::distributePointData
void distributePointData(List< T > &lst) const
Distribute list of point data.
Definition: mapDistributePolyMesh.H:236
zeroGradientFvPatchFields.H
Foam::argList::optionReadIfPresent
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::fvMeshDistribute
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Definition: fvMeshDistribute.H:70
Foam::rmDir
bool rmDir(const fileName &)
Remove a dirctory and its contents.
Definition: POSIX.C:974
Foam::argList::globalCaseName
const fileName & globalCaseName() const
Return case name.
Definition: argListI.H:48
Foam::readFields
void readFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool readOldTime)
Definition: readFields.C:33
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::parFvFieldReconstructor
Finite volume reconstructor for volume and surface fields.
Definition: parFvFieldReconstructor.H:57
Foam::regIOobject::timeStampMaster
@ timeStampMaster
Definition: regIOobject.H:71