decomposePar.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-2017 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  decomposePar
26 
27 Description
28  Automatically decomposes a mesh and fields of a case for parallel
29  execution of OpenFOAM.
30 
31 Usage
32  \b decomposePar [OPTION]
33 
34  Options:
35  - \par -cellDist
36  Write the cell distribution as a labelList, for use with 'manual'
37  decomposition method or as a volScalarField for post-processing.
38 
39  - \par -region <regionName> \n
40  Decompose named region. Does not check for existence of processor*.
41 
42  - \par -allRegions \n
43  Decompose all regions in regionProperties. Does not check for
44  existence of processor*.
45 
46  - \par -copyZero \n
47  Copy \a 0 directory to processor* rather than decompose the fields.
48 
49  - \par -copyUniform \n
50  Copy any \a uniform directories too.
51 
52  - \par -constant
53 
54  - \par -time xxx:yyy \n
55  Override controlDict settings and decompose selected times. Does not
56  re-decompose the mesh i.e. does not handle moving mesh or changing
57  mesh cases.
58 
59  - \par -fields \n
60  Use existing geometry decomposition and convert fields only.
61 
62  - \par -noSets \n
63  Skip decomposing cellSets, faceSets, pointSets.
64 
65  - \par -force \n
66  Remove any existing \a processor subdirectories before decomposing the
67  geometry.
68 
69  - \par -ifRequired \n
70  Only decompose the geometry if the number of domains has changed from a
71  previous decomposition. No \a processor subdirectories will be removed
72  unless the \a -force option is also specified. This option can be used
73  to avoid redundant geometry decomposition (eg, in scripts), but should
74  be used with caution when the underlying (serial) geometry or the
75  decomposition method etc. have been changed between decompositions.
76 
77  - \par -dict <filename>
78  Specify alternative dictionary for the decomposition.
79 
80 \*---------------------------------------------------------------------------*/
81 
82 #include "OSspecific.H"
83 #include "fvCFD.H"
84 #include "IOobjectList.H"
85 #include "domainDecomposition.H"
86 #include "labelIOField.H"
87 #include "labelFieldIOField.H"
88 #include "scalarIOField.H"
89 #include "scalarFieldIOField.H"
90 #include "vectorIOField.H"
91 #include "vectorFieldIOField.H"
92 #include "sphericalTensorIOField.H"
94 #include "symmTensorIOField.H"
95 #include "symmTensorFieldIOField.H"
96 #include "tensorIOField.H"
97 #include "tensorFieldIOField.H"
98 #include "pointFields.H"
99 #include "regionProperties.H"
100 
101 #include "readFields.H"
102 #include "dimFieldDecomposer.H"
103 #include "fvFieldDecomposer.H"
104 #include "pointFieldDecomposer.H"
106 #include "decompositionModel.H"
107 #include "collatedFileOperation.H"
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110 
111 namespace Foam
112 {
113 
114 const labelIOList& procAddressing
115 (
116  const PtrList<fvMesh>& procMeshList,
117  const label proci,
118  const word& name,
119  PtrList<labelIOList>& procAddressingList
120 )
121 {
122  const fvMesh& procMesh = procMeshList[proci];
123 
124  if (!procAddressingList.set(proci))
125  {
126  procAddressingList.set
127  (
128  proci,
129  new labelIOList
130  (
131  IOobject
132  (
133  name,
134  procMesh.facesInstance(),
135  procMesh.meshSubDir,
136  procMesh,
139  false
140  )
141  )
142  );
143  }
144  return procAddressingList[proci];
145 }
146 
147 
148 void decomposeUniform
149 (
150  const bool copyUniform,
151  const domainDecomposition& mesh,
152  const Time& processorDb,
153  const word& regionDir = word::null
154 )
155 {
156  const Time& runTime = mesh.time();
157 
158  // Any uniform data to copy/link?
159  const fileName uniformDir(regionDir/"uniform");
160 
161  if (fileHandler().isDir(runTime.timePath()/uniformDir))
162  {
163  Info<< "Detected additional non-decomposed files in "
164  << runTime.timePath()/uniformDir
165  << endl;
166 
167  const fileName timePath =
168  fileHandler().filePath(processorDb.timePath());
169 
170  if (copyUniform || mesh.distributed())
171  {
172  if (!fileHandler().exists(timePath/uniformDir))
173  {
174  fileHandler().cp
175  (
176  runTime.timePath()/uniformDir,
177  timePath/uniformDir
178  );
179  }
180  }
181  else
182  {
183  // link with relative paths
184  string parentPath = string("..")/"..";
185 
186  if (regionDir != word::null)
187  {
188  parentPath = parentPath/"..";
189  }
190 
191  fileName currentDir(cwd());
192  chDir(timePath);
193 
194  if (!fileHandler().exists(uniformDir))
195  {
196  fileHandler().ln
197  (
198  parentPath/runTime.timeName()/uniformDir,
199  uniformDir
200  );
201  }
202  chDir(currentDir);
203  }
204  }
205 }
206 
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 int main(int argc, char *argv[])
213 {
215  (
216  "decompose a mesh and fields of a case for parallel execution"
217  );
218 
220  #include "addRegionOption.H"
222  (
223  "allRegions",
224  "operate on all regions in regionProperties"
225  );
227  (
228  "cellDist",
229  "write cell distribution as a labelList - for use with 'manual' "
230  "decomposition method or as a volScalarField for post-processing."
231  );
233  (
234  "copyZero",
235  "Copy \a 0 directory to processor* rather than decompose the fields"
236  );
238  (
239  "copyUniform",
240  "copy any uniform/ directories too"
241  );
243  (
244  "fields",
245  "use existing geometry decomposition and convert fields only"
246  );
248  (
249  "noSets",
250  "skip decomposing cellSets, faceSets, pointSets"
251  );
253  (
254  "force",
255  "remove existing processor*/ subdirs before decomposing the geometry"
256  );
258  (
259  "ifRequired",
260  "only decompose geometry if the number of domains has changed"
261  );
262 
264  (
265  "dict",
266  "dictionary file name",
267  "specify alternative decomposition dictionary"
268  );
269 
270  // Include explicit constant options, have zero from time range
271  timeSelector::addOptions(true, false);
272 
273  #include "setRootCase.H"
274 
275  bool allRegions = args.optionFound("allRegions");
276  bool writeCellDist = args.optionFound("cellDist");
277  bool copyZero = args.optionFound("copyZero");
278  bool copyUniform = args.optionFound("copyUniform");
279  bool decomposeFieldsOnly = args.optionFound("fields");
280  bool decomposeSets = !args.optionFound("noSets");
281  bool forceOverwrite = args.optionFound("force");
282  bool ifRequiredDecomposition = args.optionFound("ifRequired");
283 
284  const word dictName("decomposeParDict");
285 
286  // Set time from database
287  #include "createTime.H"
288 
289  fileName dictPath;
290 
291  // Check if the dictionary is specified on the command-line
292  if (args.optionFound("dict"))
293  {
294  dictPath = args["dict"];
295 
296  dictPath =
297  (
298  isDir(dictPath)
299  ? dictPath/dictName
300  : dictPath
301  );
302  }
303  else
304  {
305  dictPath = runTime.path()/"system"/dictName;
306  }
307 
308  // Allow override of time
310 
311  wordList regionNames;
312  wordList regionDirs;
313  if (allRegions)
314  {
315  Info<< "Decomposing all regions in regionProperties" << nl << endl;
316  regionProperties rp(runTime);
317  forAllConstIter(HashTable<wordList>, rp, iter)
318  {
319  const wordList& regions = iter();
320  forAll(regions, i)
321  {
322  if (findIndex(regionNames, regions[i]) == -1)
323  {
324  regionNames.append(regions[i]);
325  }
326  }
327  }
328  regionDirs = regionNames;
329  }
330  else
331  {
332  word regionName;
333  if (args.optionReadIfPresent("region", regionName))
334  {
335  regionNames = wordList(1, regionName);
336  regionDirs = regionNames;
337  }
338  else
339  {
340  regionNames = wordList(1, fvMesh::defaultRegion);
341  regionDirs = wordList(1, word::null);
342  }
343  }
344 
345 
346 
347  forAll(regionNames, regioni)
348  {
349  const word& regionName = regionNames[regioni];
350  const word& regionDir = regionDirs[regioni];
351 
352  Info<< "\n\nDecomposing mesh " << regionName << nl << endl;
353 
354 
355  // Determine the existing processor count directly
356  label nProcs = fileHandler().nProcs(runTime.path(), regionDir);
357 
358  // Get requested numberOfSubdomains. Note: have no mesh yet so
359  // cannot use decompositionModel::New
360  const label nDomains = readLabel
361  (
362  IOdictionary
363  (
364  IOobject
365  (
366  dictName,
367  runTime.time().system(),
368  regionDir, // use region if non-standard
369  runTime,
372  false
373  )
374  ).lookup("numberOfSubdomains")
375  );
376 
377  if (decomposeFieldsOnly)
378  {
379  // Sanity check on previously decomposed case
380  if (nProcs != nDomains)
381  {
383  << "Specified -fields, but the case was decomposed with "
384  << nProcs << " domains"
385  << nl
386  << "instead of " << nDomains
387  << " domains as specified in " << dictName
388  << nl
389  << exit(FatalError);
390  }
391  }
392  else if (nProcs)
393  {
394  bool procDirsProblem = true;
395 
396  if (ifRequiredDecomposition && nProcs == nDomains)
397  {
398  // we can reuse the decomposition
399  decomposeFieldsOnly = true;
400  procDirsProblem = false;
401  forceOverwrite = false;
402 
403  Info<< "Using existing processor directories" << nl;
404  }
405 
406  if (forceOverwrite)
407  {
408  Info<< "Removing " << nProcs
409  << " existing processor directories" << endl;
410 
411  // Remove existing processors directory
412  const fileName procDir(runTime.path()/word("processors"));
413  if (fileHandler().exists(procDir))
414  {
415  fileHandler().rmDir(procDir);
416  }
417 
418  // Remove existing processor directories
419  // reverse order to avoid gaps if someone interrupts the process
420  for (label proci = nProcs-1; proci >= 0; --proci)
421  {
422  const fileName procDir
423  (
424  runTime.path()/(word("processor") + name(proci))
425  );
426 
427  if (fileHandler().exists(procDir))
428  {
429  fileHandler().rmDir(procDir);
430  }
431  }
432 
433  procDirsProblem = false;
434  }
435 
436  if (procDirsProblem)
437  {
439  << "Case is already decomposed with " << nProcs
440  << " domains, use the -force option or manually" << nl
441  << "remove processor directories before decomposing. e.g.,"
442  << nl
443  << " rm -rf " << runTime.path().c_str() << "/processor*"
444  << nl
445  << exit(FatalError);
446  }
447  }
448 
449  Info<< "Create mesh" << endl;
450  domainDecomposition mesh
451  (
452  IOobject
453  (
454  regionName,
455  runTime.timeName(),
456  runTime,
459  false
460  ),
461  dictPath
462  );
463 
464  // Decompose the mesh
465  if (!decomposeFieldsOnly)
466  {
467  // Disable buffering when writing mesh since we need to read
468  // it later on when decomposing the fields
469  float bufSz =
472 
473  mesh.decomposeMesh(dictPath);
474 
475  mesh.writeDecomposition(decomposeSets);
476 
477  if (writeCellDist)
478  {
479  const labelList& procIds = mesh.cellToProc();
480 
481  // Write the decomposition as labelList for use with 'manual'
482  // decomposition method.
483  labelIOList cellDecomposition
484  (
485  IOobject
486  (
487  "cellDecomposition",
488  mesh.facesInstance(),
489  mesh,
492  false
493  ),
494  procIds
495  );
496  cellDecomposition.write();
497 
498  Info<< nl << "Wrote decomposition to "
499  << cellDecomposition.objectPath()
500  << " for use in manual decomposition." << endl;
501 
502  // Write as volScalarField for postprocessing.
503  volScalarField cellDist
504  (
505  IOobject
506  (
507  "cellDist",
508  runTime.timeName(),
509  mesh,
512  ),
513  mesh,
514  dimensionedScalar("cellDist", dimless, 0)
515  );
516 
517  forAll(procIds, celli)
518  {
519  cellDist[celli] = procIds[celli];
520  }
521 
522  cellDist.write();
523 
524  Info<< nl << "Wrote decomposition as volScalarField to "
525  << cellDist.name() << " for use in postprocessing."
526  << endl;
527  }
528 
530  bufSz;
531  }
532 
533 
534  if (copyZero)
535  {
536  // Copy the 0 directory into each of the processor directories
537  fileName prevTimePath;
538  for (label proci = 0; proci < mesh.nProcs(); proci++)
539  {
540  Time processorDb
541  (
543  args.rootPath(),
544  args.caseName()/fileName(word("processor") + name(proci))
545  );
546  processorDb.setTime(runTime);
547 
548  if (fileHandler().isDir(runTime.timePath()))
549  {
550  // Get corresponding directory name (to handle processors/)
551  const fileName timePath
552  (
553  fileHandler().objectPath
554  (
555  IOobject
556  (
557  "",
558  processorDb.timeName(),
559  processorDb
560  ),
561  word::null
562  )
563  );
564 
565  if (timePath != prevTimePath)
566  {
567  Info<< "Processor " << proci
568  << ": copying " << runTime.timePath() << nl
569  << " to " << timePath << endl;
570  fileHandler().cp(runTime.timePath(), timePath);
571 
572  prevTimePath = timePath;
573  }
574  }
575  }
576  }
577  else
578  {
579  // Decompose the field files
580 
581  // Cached processor meshes and maps. These are only preserved if
582  // running with multiple times.
583  PtrList<Time> processorDbList(mesh.nProcs());
584  PtrList<fvMesh> procMeshList(mesh.nProcs());
585  PtrList<labelIOList> faceProcAddressingList(mesh.nProcs());
586  PtrList<labelIOList> cellProcAddressingList(mesh.nProcs());
587  PtrList<labelIOList> boundaryProcAddressingList(mesh.nProcs());
588  PtrList<fvFieldDecomposer> fieldDecomposerList(mesh.nProcs());
589  PtrList<dimFieldDecomposer> dimFieldDecomposerList(mesh.nProcs());
590  PtrList<labelIOList> pointProcAddressingList(mesh.nProcs());
591  PtrList<pointFieldDecomposer> pointFieldDecomposerList
592  (
593  mesh.nProcs()
594  );
595 
596 
597  // Loop over all times
598  forAll(times, timeI)
599  {
600  runTime.setTime(times[timeI], timeI);
601 
602  Info<< "Time = " << runTime.timeName() << endl;
603 
604  // Search for list of objects for this time
605  IOobjectList objects(mesh, runTime.timeName());
606 
607 
608  // Construct the vol fields
609  // ~~~~~~~~~~~~~~~~~~~~~~~~
610  PtrList<volScalarField> volScalarFields;
611  readFields(mesh, objects, volScalarFields);
612  PtrList<volVectorField> volVectorFields;
613  readFields(mesh, objects, volVectorFields);
614  PtrList<volSphericalTensorField> volSphericalTensorFields;
615  readFields(mesh, objects, volSphericalTensorFields);
616  PtrList<volSymmTensorField> volSymmTensorFields;
617  readFields(mesh, objects, volSymmTensorFields);
618  PtrList<volTensorField> volTensorFields;
619  readFields(mesh, objects, volTensorFields);
620 
621 
622  // Construct the dimensioned fields
623  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
624  PtrList<DimensionedField<scalar, volMesh>> dimScalarFields;
625  readFields(mesh, objects, dimScalarFields);
626  PtrList<DimensionedField<vector, volMesh>> dimVectorFields;
627  readFields(mesh, objects, dimVectorFields);
628  PtrList<DimensionedField<sphericalTensor, volMesh>>
629  dimSphericalTensorFields;
630  readFields(mesh, objects, dimSphericalTensorFields);
631  PtrList<DimensionedField<symmTensor, volMesh>>
632  dimSymmTensorFields;
633  readFields(mesh, objects, dimSymmTensorFields);
634  PtrList<DimensionedField<tensor, volMesh>> dimTensorFields;
635  readFields(mesh, objects, dimTensorFields);
636 
637 
638  // Construct the surface fields
639  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
640  PtrList<surfaceScalarField> surfaceScalarFields;
641  readFields(mesh, objects, surfaceScalarFields);
642  PtrList<surfaceVectorField> surfaceVectorFields;
643  readFields(mesh, objects, surfaceVectorFields);
644  PtrList<surfaceSphericalTensorField>
645  surfaceSphericalTensorFields;
646  readFields(mesh, objects, surfaceSphericalTensorFields);
647  PtrList<surfaceSymmTensorField> surfaceSymmTensorFields;
648  readFields(mesh, objects, surfaceSymmTensorFields);
649  PtrList<surfaceTensorField> surfaceTensorFields;
650  readFields(mesh, objects, surfaceTensorFields);
651 
652 
653  // Construct the point fields
654  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
655  const pointMesh& pMesh = pointMesh::New(mesh);
656 
657  PtrList<pointScalarField> pointScalarFields;
658  readFields(pMesh, objects, pointScalarFields);
659  PtrList<pointVectorField> pointVectorFields;
660  readFields(pMesh, objects, pointVectorFields);
661  PtrList<pointSphericalTensorField> pointSphericalTensorFields;
662  readFields(pMesh, objects, pointSphericalTensorFields);
663  PtrList<pointSymmTensorField> pointSymmTensorFields;
664  readFields(pMesh, objects, pointSymmTensorFields);
665  PtrList<pointTensorField> pointTensorFields;
666  readFields(pMesh, objects, pointTensorFields);
667 
668 
669  // Construct the Lagrangian fields
670  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
671 
672  fileNameList cloudDirs
673  (
675  (
676  runTime.timePath()/cloud::prefix,
678  )
679  );
680 
681  // Particles
682  PtrList<Cloud<indexedParticle>> lagrangianPositions
683  (
684  cloudDirs.size()
685  );
686  // Particles per cell
687  PtrList<List<SLList<indexedParticle*>*>> cellParticles
688  (
689  cloudDirs.size()
690  );
691 
692  PtrList<PtrList<labelIOField>> lagrangianLabelFields
693  (
694  cloudDirs.size()
695  );
696  PtrList<PtrList<labelFieldCompactIOField>>
697  lagrangianLabelFieldFields
698  (
699  cloudDirs.size()
700  );
701  PtrList<PtrList<scalarIOField>> lagrangianScalarFields
702  (
703  cloudDirs.size()
704  );
705  PtrList<PtrList<scalarFieldCompactIOField>>
706  lagrangianScalarFieldFields
707  (
708  cloudDirs.size()
709  );
710  PtrList<PtrList<vectorIOField>> lagrangianVectorFields
711  (
712  cloudDirs.size()
713  );
714  PtrList<PtrList<vectorFieldCompactIOField>>
715  lagrangianVectorFieldFields
716  (
717  cloudDirs.size()
718  );
719  PtrList<PtrList<sphericalTensorIOField>>
720  lagrangianSphericalTensorFields
721  (
722  cloudDirs.size()
723  );
724  PtrList<PtrList<sphericalTensorFieldCompactIOField>>
725  lagrangianSphericalTensorFieldFields(cloudDirs.size());
726  PtrList<PtrList<symmTensorIOField>> lagrangianSymmTensorFields
727  (
728  cloudDirs.size()
729  );
730  PtrList<PtrList<symmTensorFieldCompactIOField>>
731  lagrangianSymmTensorFieldFields
732  (
733  cloudDirs.size()
734  );
735  PtrList<PtrList<tensorIOField>> lagrangianTensorFields
736  (
737  cloudDirs.size()
738  );
739  PtrList<PtrList<tensorFieldCompactIOField>>
740  lagrangianTensorFieldFields
741  (
742  cloudDirs.size()
743  );
744 
745  label cloudI = 0;
746 
747  forAll(cloudDirs, i)
748  {
749  IOobjectList sprayObjs
750  (
751  mesh,
752  runTime.timeName(),
753  cloud::prefix/cloudDirs[i],
756  false
757  );
758 
759  IOobject* positionsPtr = sprayObjs.lookup
760  (
761  word("positions")
762  );
763 
764  if (positionsPtr)
765  {
766  // Read lagrangian particles
767  // ~~~~~~~~~~~~~~~~~~~~~~~~~
768 
769  Info<< "Identified lagrangian data set: "
770  << cloudDirs[i] << endl;
771 
772  lagrangianPositions.set
773  (
774  cloudI,
775  new Cloud<indexedParticle>
776  (
777  mesh,
778  cloudDirs[i],
779  false
780  )
781  );
782 
783 
784  // Sort particles per cell
785  // ~~~~~~~~~~~~~~~~~~~~~~~
786 
787  cellParticles.set
788  (
789  cloudI,
790  new List<SLList<indexedParticle*>*>
791  (
792  mesh.nCells(),
793  static_cast<SLList<indexedParticle*>*>(nullptr)
794  )
795  );
796 
797  label i = 0;
798 
799  forAllIter
800  (
801  Cloud<indexedParticle>,
802  lagrangianPositions[cloudI],
803  iter
804  )
805  {
806  iter().index() = i++;
807 
808  label celli = iter().cell();
809 
810  // Check
811  if (celli < 0 || celli >= mesh.nCells())
812  {
814  << "Illegal cell number " << celli
815  << " for particle with index "
816  << iter().index()
817  << " at position "
818  << iter().position() << nl
819  << "Cell number should be between 0 and "
820  << mesh.nCells()-1 << nl
821  << "On this mesh the particle should"
822  << " be in cell "
823  << mesh.findCell(iter().position())
824  << exit(FatalError);
825  }
826 
827  if (!cellParticles[cloudI][celli])
828  {
829  cellParticles[cloudI][celli] =
830  new SLList<indexedParticle*>();
831  }
832 
833  cellParticles[cloudI][celli]->append(&iter());
834  }
835 
836  // Read fields
837  // ~~~~~~~~~~~
838 
839  IOobjectList lagrangianObjects
840  (
841  mesh,
842  runTime.timeName(),
843  cloud::prefix/cloudDirs[cloudI],
846  false
847  );
848 
850  (
851  cloudI,
852  lagrangianObjects,
853  lagrangianLabelFields
854  );
855 
857  (
858  cloudI,
859  lagrangianObjects,
860  lagrangianLabelFieldFields
861  );
862 
864  (
865  cloudI,
866  lagrangianObjects,
867  lagrangianScalarFields
868  );
869 
871  (
872  cloudI,
873  lagrangianObjects,
874  lagrangianScalarFieldFields
875  );
876 
878  (
879  cloudI,
880  lagrangianObjects,
881  lagrangianVectorFields
882  );
883 
885  (
886  cloudI,
887  lagrangianObjects,
888  lagrangianVectorFieldFields
889  );
890 
892  (
893  cloudI,
894  lagrangianObjects,
895  lagrangianSphericalTensorFields
896  );
897 
899  (
900  cloudI,
901  lagrangianObjects,
902  lagrangianSphericalTensorFieldFields
903  );
904 
906  (
907  cloudI,
908  lagrangianObjects,
909  lagrangianSymmTensorFields
910  );
911 
913  (
914  cloudI,
915  lagrangianObjects,
916  lagrangianSymmTensorFieldFields
917  );
918 
920  (
921  cloudI,
922  lagrangianObjects,
923  lagrangianTensorFields
924  );
925 
927  (
928  cloudI,
929  lagrangianObjects,
930  lagrangianTensorFieldFields
931  );
932 
933  cloudI++;
934  }
935  }
936 
937  lagrangianPositions.setSize(cloudI);
938  cellParticles.setSize(cloudI);
939  lagrangianLabelFields.setSize(cloudI);
940  lagrangianLabelFieldFields.setSize(cloudI);
941  lagrangianScalarFields.setSize(cloudI);
942  lagrangianScalarFieldFields.setSize(cloudI);
943  lagrangianVectorFields.setSize(cloudI);
944  lagrangianVectorFieldFields.setSize(cloudI);
945  lagrangianSphericalTensorFields.setSize(cloudI);
946  lagrangianSphericalTensorFieldFields.setSize(cloudI);
947  lagrangianSymmTensorFields.setSize(cloudI);
948  lagrangianSymmTensorFieldFields.setSize(cloudI);
949  lagrangianTensorFields.setSize(cloudI);
950  lagrangianTensorFieldFields.setSize(cloudI);
951 
952  Info<< endl;
953 
954  // split the fields over processors
955  for (label proci = 0; proci < mesh.nProcs(); proci++)
956  {
957  Info<< "Processor " << proci << ": field transfer" << endl;
958 
959 
960  // open the database
961  if (!processorDbList.set(proci))
962  {
963  processorDbList.set
964  (
965  proci,
966  new Time
967  (
969  args.rootPath(),
970  args.caseName()
971  /fileName(word("processor") + name(proci))
972  )
973  );
974  }
975  Time& processorDb = processorDbList[proci];
976 
977 
978  processorDb.setTime(runTime);
979 
980  // read the mesh
981  if (!procMeshList.set(proci))
982  {
983  procMeshList.set
984  (
985  proci,
986  new fvMesh
987  (
988  IOobject
989  (
990  regionName,
991  processorDb.timeName(),
992  processorDb
993  )
994  )
995  );
996  }
997  const fvMesh& procMesh = procMeshList[proci];
998 
999  const labelIOList& faceProcAddressing = procAddressing
1000  (
1001  procMeshList,
1002  proci,
1003  "faceProcAddressing",
1004  faceProcAddressingList
1005  );
1006 
1007  const labelIOList& cellProcAddressing = procAddressing
1008  (
1009  procMeshList,
1010  proci,
1011  "cellProcAddressing",
1012  cellProcAddressingList
1013  );
1014 
1015  const labelIOList& boundaryProcAddressing = procAddressing
1016  (
1017  procMeshList,
1018  proci,
1019  "boundaryProcAddressing",
1020  boundaryProcAddressingList
1021  );
1022 
1023 
1024  // FV fields
1025  {
1026  if (!fieldDecomposerList.set(proci))
1027  {
1028  fieldDecomposerList.set
1029  (
1030  proci,
1031  new fvFieldDecomposer
1032  (
1033  mesh,
1034  procMesh,
1035  faceProcAddressing,
1036  cellProcAddressing,
1037  boundaryProcAddressing
1038  )
1039  );
1040  }
1041  const fvFieldDecomposer& fieldDecomposer =
1042  fieldDecomposerList[proci];
1043 
1044  fieldDecomposer.decomposeFields(volScalarFields);
1045  fieldDecomposer.decomposeFields(volVectorFields);
1046  fieldDecomposer.decomposeFields
1047  (
1048  volSphericalTensorFields
1049  );
1050  fieldDecomposer.decomposeFields(volSymmTensorFields);
1051  fieldDecomposer.decomposeFields(volTensorFields);
1052 
1053  fieldDecomposer.decomposeFields(surfaceScalarFields);
1054  fieldDecomposer.decomposeFields(surfaceVectorFields);
1055  fieldDecomposer.decomposeFields
1056  (
1057  surfaceSphericalTensorFields
1058  );
1059  fieldDecomposer.decomposeFields
1060  (
1061  surfaceSymmTensorFields
1062  );
1063  fieldDecomposer.decomposeFields(surfaceTensorFields);
1064 
1065  if (times.size() == 1)
1066  {
1067  // Clear cached decomposer
1068  fieldDecomposerList.set(proci, nullptr);
1069  }
1070  }
1071 
1072  // Dimensioned fields
1073  {
1074  if (!dimFieldDecomposerList.set(proci))
1075  {
1076  dimFieldDecomposerList.set
1077  (
1078  proci,
1079  new dimFieldDecomposer
1080  (
1081  mesh,
1082  procMesh,
1083  faceProcAddressing,
1084  cellProcAddressing
1085  )
1086  );
1087  }
1088  const dimFieldDecomposer& dimDecomposer =
1089  dimFieldDecomposerList[proci];
1090 
1091  dimDecomposer.decomposeFields(dimScalarFields);
1092  dimDecomposer.decomposeFields(dimVectorFields);
1093  dimDecomposer.decomposeFields(dimSphericalTensorFields);
1094  dimDecomposer.decomposeFields(dimSymmTensorFields);
1095  dimDecomposer.decomposeFields(dimTensorFields);
1096 
1097  if (times.size() == 1)
1098  {
1099  dimFieldDecomposerList.set(proci, nullptr);
1100  }
1101  }
1102 
1103 
1104  // Point fields
1105  if
1106  (
1107  pointScalarFields.size()
1108  || pointVectorFields.size()
1109  || pointSphericalTensorFields.size()
1110  || pointSymmTensorFields.size()
1111  || pointTensorFields.size()
1112  )
1113  {
1114  const labelIOList& pointProcAddressing = procAddressing
1115  (
1116  procMeshList,
1117  proci,
1118  "pointProcAddressing",
1119  pointProcAddressingList
1120  );
1121 
1122  const pointMesh& procPMesh = pointMesh::New(procMesh);
1123 
1124  if (!pointFieldDecomposerList.set(proci))
1125  {
1126  pointFieldDecomposerList.set
1127  (
1128  proci,
1129  new pointFieldDecomposer
1130  (
1131  pMesh,
1132  procPMesh,
1133  pointProcAddressing,
1134  boundaryProcAddressing
1135  )
1136  );
1137  }
1138  const pointFieldDecomposer& pointDecomposer =
1139  pointFieldDecomposerList[proci];
1140 
1141  pointDecomposer.decomposeFields(pointScalarFields);
1142  pointDecomposer.decomposeFields(pointVectorFields);
1143  pointDecomposer.decomposeFields
1144  (
1145  pointSphericalTensorFields
1146  );
1147  pointDecomposer.decomposeFields(pointSymmTensorFields);
1148  pointDecomposer.decomposeFields(pointTensorFields);
1149 
1150 
1151  if (times.size() == 1)
1152  {
1153  pointProcAddressingList.set(proci, nullptr);
1154  pointFieldDecomposerList.set(proci, nullptr);
1155  }
1156  }
1157 
1158 
1159  // If there is lagrangian data write it out
1160  forAll(lagrangianPositions, cloudI)
1161  {
1162  if (lagrangianPositions[cloudI].size())
1163  {
1164  lagrangianFieldDecomposer fieldDecomposer
1165  (
1166  mesh,
1167  procMesh,
1168  faceProcAddressing,
1169  cellProcAddressing,
1170  cloudDirs[cloudI],
1171  lagrangianPositions[cloudI],
1172  cellParticles[cloudI]
1173  );
1174 
1175  // Lagrangian fields
1176  {
1177  fieldDecomposer.decomposeFields
1178  (
1179  cloudDirs[cloudI],
1180  lagrangianLabelFields[cloudI]
1181  );
1182  fieldDecomposer.decomposeFieldFields
1183  (
1184  cloudDirs[cloudI],
1185  lagrangianLabelFieldFields[cloudI]
1186  );
1187  fieldDecomposer.decomposeFields
1188  (
1189  cloudDirs[cloudI],
1190  lagrangianScalarFields[cloudI]
1191  );
1192  fieldDecomposer.decomposeFieldFields
1193  (
1194  cloudDirs[cloudI],
1195  lagrangianScalarFieldFields[cloudI]
1196  );
1197  fieldDecomposer.decomposeFields
1198  (
1199  cloudDirs[cloudI],
1200  lagrangianVectorFields[cloudI]
1201  );
1202  fieldDecomposer.decomposeFieldFields
1203  (
1204  cloudDirs[cloudI],
1205  lagrangianVectorFieldFields[cloudI]
1206  );
1207  fieldDecomposer.decomposeFields
1208  (
1209  cloudDirs[cloudI],
1210  lagrangianSphericalTensorFields[cloudI]
1211  );
1212  fieldDecomposer.decomposeFieldFields
1213  (
1214  cloudDirs[cloudI],
1215  lagrangianSphericalTensorFieldFields[cloudI]
1216  );
1217  fieldDecomposer.decomposeFields
1218  (
1219  cloudDirs[cloudI],
1220  lagrangianSymmTensorFields[cloudI]
1221  );
1222  fieldDecomposer.decomposeFieldFields
1223  (
1224  cloudDirs[cloudI],
1225  lagrangianSymmTensorFieldFields[cloudI]
1226  );
1227  fieldDecomposer.decomposeFields
1228  (
1229  cloudDirs[cloudI],
1230  lagrangianTensorFields[cloudI]
1231  );
1232  fieldDecomposer.decomposeFieldFields
1233  (
1234  cloudDirs[cloudI],
1235  lagrangianTensorFieldFields[cloudI]
1236  );
1237  }
1238  }
1239  }
1240 
1241  // Decompose the "uniform" directory in the time region
1242  // directory
1243  decomposeUniform(copyUniform, mesh, processorDb, regionDir);
1244 
1245  // For the first region of a multi-region case additionally
1246  // decompose the "uniform" directory in the time directory
1247  if (regionNames.size() > 1 && regioni == 0)
1248  {
1249  decomposeUniform(copyUniform, mesh, processorDb);
1250  }
1251 
1252  // We have cached all the constant mesh data for the current
1253  // processor. This is only important if running with
1254  // multiple times, otherwise it is just extra storage.
1255  if (times.size() == 1)
1256  {
1257  boundaryProcAddressingList.set(proci, nullptr);
1258  cellProcAddressingList.set(proci, nullptr);
1259  faceProcAddressingList.set(proci, nullptr);
1260  procMeshList.set(proci, nullptr);
1261  processorDbList.set(proci, nullptr);
1262  }
1263  }
1264  }
1265  }
1266  }
1267 
1268  Info<< "\nEnd\n" << endl;
1269 
1270  return 0;
1271 }
1272 
1273 
1274 // ************************************************************************* //
List< instant > instantList
List of instants.
Definition: instantList.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual bool rmDir(const fileName &) const =0
Remove a dirctory and its contents.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const fileName & rootPath() const
Return root path.
Definition: argListI.H:36
Foam::word regionName
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
Definition: POSIX.C:291
static void noParallel()
Remove the parallel options.
Definition: argList.C:148
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
virtual bool exists(const fileName &, const bool checkGzip=true, const bool followLink=true) const =0
Does the name exist (as DIRECTORY or FILE) in the file system?
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
static void readFieldFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< CompactIOField< Field< Type >, Type >> > &lagrangianFields)
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
dynamicFvMesh & mesh
virtual bool ln(const fileName &src, const fileName &dst) const =0
Create a softlink. dst should not exist. Returns true if.
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:536
fileName dictPath
regionProperties rp(runTime)
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:95
static const word null
An empty word.
Definition: word.H:77
List< label > labelList
A List of labels.
Definition: labelList.H:56
const word dictName("particleTrackDict")
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times (as select0)
Definition: timeSelector.C:288
virtual fileName filePath(const bool checkGlobal, const IOobject &, const word &typeName) const =0
Search for an object. checkGlobal : also check undecomposed case.
const fileOperation & fileHandler()
Get current file handler.
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:201
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:42
static const char nl
Definition: Ostream.H:262
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static float maxThreadFileBufferSize
Max size of thread buffer size. This is the overall size of.
fileName path() const
Return the path to the caseName.
Definition: argListI.H:60
List< word > wordList
A List of words.
Definition: fileName.H:54
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:337
fileName cwd()
Return current working directory path name.
Definition: POSIX.C:248
fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true, const bool followLink=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:641
messageStream Info
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:85
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:126
static void readFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< IOField< Type >>> &lagrangianFields)
Foam::argList args(argc, argv)
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
bool exists(const fileName &, const bool checkGzip=true, const bool followLink=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:517
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Namespace for OpenFOAM.