decomposePar.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 namespace Foam
111 {
112 
113 const labelIOList& procAddressing
114 (
115  const PtrList<fvMesh>& procMeshList,
116  const label proci,
117  const word& name,
118  PtrList<labelIOList>& procAddressingList
119 )
120 {
121  const fvMesh& procMesh = procMeshList[proci];
122 
123  if (!procAddressingList.set(proci))
124  {
125  procAddressingList.set
126  (
127  proci,
128  new labelIOList
129  (
130  IOobject
131  (
132  name,
133  procMesh.facesInstance(),
134  procMesh.meshSubDir,
135  procMesh,
138  false
139  )
140  )
141  );
142  }
143  return procAddressingList[proci];
144 }
145 
146 
147 void decomposeUniform
148 (
149  const bool copyUniform,
150  const domainDecomposition& mesh,
151  const Time& processorDb,
152  const word& regionDir = word::null
153 )
154 {
155  const Time& runTime = mesh.time();
156 
157  // Any uniform data to copy/link?
158  const fileName uniformDir(regionDir/"uniform");
159 
160  if (fileHandler().isDir(runTime.timePath()/uniformDir))
161  {
162  Info<< "Detected additional non-decomposed files in "
163  << runTime.timePath()/uniformDir
164  << endl;
165 
166  const fileName timePath =
167  fileHandler().filePath(processorDb.timePath());
168 
169  if (copyUniform || mesh.distributed())
170  {
171  if (!fileHandler().exists(timePath/uniformDir))
172  {
173  fileHandler().cp
174  (
175  runTime.timePath()/uniformDir,
176  timePath/uniformDir
177  );
178  }
179  }
180  else
181  {
182  // link with relative paths
183  string parentPath = string("..")/"..";
184 
185  if (regionDir != word::null)
186  {
187  parentPath = parentPath/"..";
188  }
189 
190  fileName currentDir(cwd());
191  chDir(timePath);
192 
193  if (!fileHandler().exists(uniformDir))
194  {
195  fileHandler().ln
196  (
197  parentPath/runTime.timeName()/uniformDir,
198  uniformDir
199  );
200  }
201  chDir(currentDir);
202  }
203  }
204 }
205 
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 int main(int argc, char *argv[])
212 {
214  (
215  "decompose a mesh and fields of a case for parallel execution"
216  );
217 
219  #include "addRegionOption.H"
220  #include "addAllRegionsOption.H"
222  (
223  "cellDist",
224  "write cell distribution as a labelList - for use with 'manual' "
225  "decomposition method or as a volScalarField for post-processing."
226  );
228  (
229  "copyZero",
230  "Copy \a 0 directory to processor* rather than decompose the fields"
231  );
233  (
234  "copyUniform",
235  "copy any uniform/ directories too"
236  );
238  (
239  "fields",
240  "use existing geometry decomposition and convert fields only"
241  );
243  (
244  "noFields",
245  "opposite of -fields; only decompose geometry"
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 region = args.optionFound("region");
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 decomposeGeomOnly = args.optionFound("noFields");
281  bool decomposeSets = !args.optionFound("noSets");
282  bool forceOverwrite = args.optionFound("force");
283  bool ifRequiredDecomposition = args.optionFound("ifRequired");
284 
285  const word dictName("decomposeParDict");
286 
287 
288  if (decomposeGeomOnly)
289  {
290  Info<< "Skipping decomposing fields"
291  << nl << endl;
292 
293  if (decomposeFieldsOnly || copyZero)
294  {
296  << "Cannot combine geometry-only decomposition (-noFields)"
297  << " with field decomposition (-noFields or -copyZero)"
298  << exit(FatalError);
299  }
300  }
301 
302 
303  // Set time from database
304  #include "createTime.H"
305 
306  // Check if the dictionary is specified on the command-line
307  fileName dictPath = fileName::null;
308  if (args.optionFound("dict"))
309  {
310  dictPath = args["dict"];
311 
312  if (!isFile(dictPath))
313  {
314  dictPath = dictPath/dictName;
315  }
316 
317  if (!isFile(dictPath))
318  {
320  << "Specified -dict " << args["dict"] << " but neither "
321  << args["dict"] << " nor " << args["dict"]/dictName
322  << " could be found" << nl << exit(FatalError);
323  }
324  }
325 
326  // Allow override of time
328 
329  const wordList regionNames(selectRegionNames(args, runTime));
330 
331  {
332  // Determine the existing processor count directly
333  label nProcs = fileHandler().nProcs(runTime.path());
334 
335  if (forceOverwrite)
336  {
337  if (region)
338  {
340  << "Cannot force the decomposition of a single region"
341  << exit(FatalError);
342  }
343 
344  Info<< "Removing " << nProcs
345  << " existing processor directories" << endl;
346 
347  // Remove existing processors directory
348  fileNameList dirs
349  (
351  (
352  runTime.path(),
354  )
355  );
356  forAllReverse(dirs, diri)
357  {
358  const fileName& d = dirs[diri];
359 
360  // Starts with 'processors'
361  if (d.find("processors") == 0)
362  {
363  if (fileHandler().exists(d))
364  {
365  fileHandler().rmDir(d);
366  }
367  }
368 
369  // Starts with 'processor'
370  if (d.find("processor") == 0)
371  {
372  // Check that integer after processor
373  fileName num(d.substr(9));
374  label proci = -1;
375  if (Foam::read(num.c_str(), proci))
376  {
377  if (fileHandler().exists(d))
378  {
379  fileHandler().rmDir(d);
380  }
381  }
382  }
383  }
384  }
385  else if (nProcs && !region && !decomposeFieldsOnly)
386  {
388  << "Case is already decomposed with " << nProcs
389  << " domains, use the -force option or manually" << nl
390  << "remove processor directories before decomposing. e.g.,"
391  << nl
392  << " rm -rf " << runTime.path().c_str() << "/processor*"
393  << nl
394  << exit(FatalError);
395  }
396  }
397 
398 
399  forAll(regionNames, regioni)
400  {
401  const word& regionName = regionNames[regioni];
402  const word& regionDir = Foam::regionDir(regionName);
403 
404  Info<< "\n\nDecomposing mesh " << regionName << nl << endl;
405 
406  // Determine the existing processor count directly
407  label nProcs = fileHandler().nProcs(runTime.path(), regionDir);
408 
409  // Get the dictionary IO
410  const IOobject dictIO
411  (
412  dictPath == fileName::null
413  ? IOobject
414  (
415  dictName,
416  runTime.time().system(),
417  regionDir, // use region if non-standard
418  runTime,
421  false
422  )
423  : IOobject
424  (
425  dictPath,
426  runTime,
427  IOobject::MUST_READ_IF_MODIFIED,
428  IOobject::NO_WRITE,
429  false
430  )
431  );
432  // Get requested numberOfSubdomains. Note: have no mesh yet so
433  // cannot use decompositionModel::New
434  const label nDomains =
435  readLabel(IOdictionary(dictIO).lookup("numberOfSubdomains"));
436 
437  // Give file handler a chance to determine the output directory
438  const_cast<fileOperation&>(fileHandler()).setNProcs(nDomains);
439 
440  if (decomposeFieldsOnly)
441  {
442  // Sanity check on previously decomposed case
443  if (nProcs != nDomains)
444  {
446  << "Specified -fields, but the case was decomposed with "
447  << nProcs << " domains"
448  << nl
449  << "instead of " << nDomains
450  << " domains as specified in " << dictName
451  << nl
452  << exit(FatalError);
453  }
454  }
455  else if (nProcs)
456  {
457  if (ifRequiredDecomposition && nProcs == nDomains)
458  {
459  // Reuse the decomposition
460  decomposeFieldsOnly = true;
461  Info<< "Using existing processor directories" << nl;
462  }
463  }
464 
465  Info<< "Create mesh" << endl;
466  domainDecomposition mesh
467  (
468  IOobject
469  (
470  regionName,
471  runTime.timeName(),
472  runTime,
475  false
476  ),
478  );
479 
480  // Decompose the mesh
481  if (!decomposeFieldsOnly)
482  {
483  mesh.decomposeMesh(dictIO.objectPath());
484 
485  mesh.writeDecomposition(decomposeSets);
486 
487  if (writeCellDist)
488  {
489  const labelList& procIds = mesh.cellToProc();
490 
491  // Write the decomposition as labelList for use with 'manual'
492  // decomposition method.
493  labelIOList cellDecomposition
494  (
495  IOobject
496  (
497  "cellDecomposition",
498  mesh.facesInstance(),
499  mesh,
502  false
503  ),
504  procIds
505  );
506  cellDecomposition.write();
507 
508  Info<< nl << "Wrote decomposition to "
509  << cellDecomposition.objectPath()
510  << " for use in manual decomposition." << endl;
511 
512  // Write as volScalarField for postprocessing.
513  volScalarField cellDist
514  (
515  IOobject
516  (
517  "cellDist",
518  runTime.timeName(),
519  mesh,
522  ),
523  mesh,
525  );
526 
527  forAll(procIds, celli)
528  {
529  cellDist[celli] = procIds[celli];
530  }
531 
532  cellDist.write();
533 
534  Info<< nl << "Wrote decomposition as volScalarField to "
535  << cellDist.name() << " for use in postprocessing."
536  << endl;
537  }
538 
539  fileHandler().flush();
540  }
541 
542 
543  if (copyZero)
544  {
545  // Copy the 0 directory into each of the processor directories
546  fileName prevTimePath;
547  for (label proci = 0; proci < mesh.nProcs(); proci++)
548  {
549  Time processorDb
550  (
552  args.rootPath(),
553  args.caseName()/fileName(word("processor") + name(proci))
554  );
555  processorDb.setTime(runTime);
556 
557  if (fileHandler().isDir(runTime.timePath()))
558  {
559  // Get corresponding directory name (to handle processors/)
560  const fileName timePath
561  (
562  fileHandler().objectPath
563  (
564  IOobject
565  (
566  "",
567  processorDb.timeName(),
568  processorDb
569  ),
570  word::null
571  )
572  );
573 
574  if (timePath != prevTimePath)
575  {
576  Info<< "Processor " << proci
577  << ": copying " << runTime.timePath() << nl
578  << " to " << timePath << endl;
579  fileHandler().cp(runTime.timePath(), timePath);
580 
581  prevTimePath = timePath;
582  }
583  }
584  }
585  }
586  else if (!decomposeGeomOnly)
587  {
588  // Decompose the field files
589 
590  // Cached processor meshes and maps. These are only preserved if
591  // running with multiple times.
592  PtrList<Time> processorDbList(mesh.nProcs());
593  PtrList<fvMesh> procMeshList(mesh.nProcs());
594  PtrList<labelIOList> faceProcAddressingList(mesh.nProcs());
595  PtrList<labelIOList> cellProcAddressingList(mesh.nProcs());
596  PtrList<labelIOList> boundaryProcAddressingList(mesh.nProcs());
597  PtrList<fvFieldDecomposer> fieldDecomposerList(mesh.nProcs());
598  PtrList<dimFieldDecomposer> dimFieldDecomposerList(mesh.nProcs());
599  PtrList<labelIOList> pointProcAddressingList(mesh.nProcs());
600  PtrList<pointFieldDecomposer> pointFieldDecomposerList
601  (
602  mesh.nProcs()
603  );
604 
605 
606  // Loop over all times
607  forAll(times, timeI)
608  {
609  runTime.setTime(times[timeI], timeI);
610 
611  Info<< "Time = " << runTime.timeName() << endl;
612 
613  // Search for list of objects for this time
614  IOobjectList objects(mesh, runTime.timeName());
615 
616 
617  // Construct the vol fields
618  // ~~~~~~~~~~~~~~~~~~~~~~~~
619  PtrList<volScalarField> volScalarFields;
620  readFields(mesh, objects, volScalarFields);
621  PtrList<volVectorField> volVectorFields;
622  readFields(mesh, objects, volVectorFields);
623  PtrList<volSphericalTensorField> volSphericalTensorFields;
624  readFields(mesh, objects, volSphericalTensorFields);
625  PtrList<volSymmTensorField> volSymmTensorFields;
626  readFields(mesh, objects, volSymmTensorFields);
627  PtrList<volTensorField> volTensorFields;
628  readFields(mesh, objects, volTensorFields);
629 
630 
631  // Construct the dimensioned fields
632  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
633  PtrList<DimensionedField<scalar, volMesh>> dimScalarFields;
634  readFields(mesh, objects, dimScalarFields);
635  PtrList<DimensionedField<vector, volMesh>> dimVectorFields;
636  readFields(mesh, objects, dimVectorFields);
637  PtrList<DimensionedField<sphericalTensor, volMesh>>
638  dimSphericalTensorFields;
639  readFields(mesh, objects, dimSphericalTensorFields);
640  PtrList<DimensionedField<symmTensor, volMesh>>
641  dimSymmTensorFields;
642  readFields(mesh, objects, dimSymmTensorFields);
643  PtrList<DimensionedField<tensor, volMesh>> dimTensorFields;
644  readFields(mesh, objects, dimTensorFields);
645 
646 
647  // Construct the surface fields
648  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
649  PtrList<surfaceScalarField> surfaceScalarFields;
650  readFields(mesh, objects, surfaceScalarFields);
651  PtrList<surfaceVectorField> surfaceVectorFields;
652  readFields(mesh, objects, surfaceVectorFields);
653  PtrList<surfaceSphericalTensorField>
654  surfaceSphericalTensorFields;
655  readFields(mesh, objects, surfaceSphericalTensorFields);
656  PtrList<surfaceSymmTensorField> surfaceSymmTensorFields;
657  readFields(mesh, objects, surfaceSymmTensorFields);
658  PtrList<surfaceTensorField> surfaceTensorFields;
659  readFields(mesh, objects, surfaceTensorFields);
660 
661 
662  // Construct the point fields
663  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
664  const pointMesh& pMesh = pointMesh::New(mesh);
665 
666  PtrList<pointScalarField> pointScalarFields;
667  readFields(pMesh, objects, pointScalarFields);
668  PtrList<pointVectorField> pointVectorFields;
669  readFields(pMesh, objects, pointVectorFields);
670  PtrList<pointSphericalTensorField> pointSphericalTensorFields;
671  readFields(pMesh, objects, pointSphericalTensorFields);
672  PtrList<pointSymmTensorField> pointSymmTensorFields;
673  readFields(pMesh, objects, pointSymmTensorFields);
674  PtrList<pointTensorField> pointTensorFields;
675  readFields(pMesh, objects, pointTensorFields);
676 
677 
678  // Construct the Lagrangian fields
679  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
680 
681  fileNameList cloudDirs
682  (
684  (
685  runTime.timePath()/cloud::prefix,
687  )
688  );
689 
690  // Particles
691  PtrList<Cloud<indexedParticle>> lagrangianPositions
692  (
693  cloudDirs.size()
694  );
695  // Particles per cell
696  PtrList<List<SLList<indexedParticle*>*>> cellParticles
697  (
698  cloudDirs.size()
699  );
700 
701  PtrList<PtrList<labelIOField>> lagrangianLabelFields
702  (
703  cloudDirs.size()
704  );
705  PtrList<PtrList<labelFieldCompactIOField>>
706  lagrangianLabelFieldFields
707  (
708  cloudDirs.size()
709  );
710  PtrList<PtrList<scalarIOField>> lagrangianScalarFields
711  (
712  cloudDirs.size()
713  );
714  PtrList<PtrList<scalarFieldCompactIOField>>
715  lagrangianScalarFieldFields
716  (
717  cloudDirs.size()
718  );
719  PtrList<PtrList<vectorIOField>> lagrangianVectorFields
720  (
721  cloudDirs.size()
722  );
723  PtrList<PtrList<vectorFieldCompactIOField>>
724  lagrangianVectorFieldFields
725  (
726  cloudDirs.size()
727  );
728  PtrList<PtrList<sphericalTensorIOField>>
729  lagrangianSphericalTensorFields
730  (
731  cloudDirs.size()
732  );
733  PtrList<PtrList<sphericalTensorFieldCompactIOField>>
734  lagrangianSphericalTensorFieldFields(cloudDirs.size());
735  PtrList<PtrList<symmTensorIOField>> lagrangianSymmTensorFields
736  (
737  cloudDirs.size()
738  );
739  PtrList<PtrList<symmTensorFieldCompactIOField>>
740  lagrangianSymmTensorFieldFields
741  (
742  cloudDirs.size()
743  );
744  PtrList<PtrList<tensorIOField>> lagrangianTensorFields
745  (
746  cloudDirs.size()
747  );
748  PtrList<PtrList<tensorFieldCompactIOField>>
749  lagrangianTensorFieldFields
750  (
751  cloudDirs.size()
752  );
753 
754  label cloudI = 0;
755 
756  forAll(cloudDirs, i)
757  {
758  IOobjectList sprayObjs
759  (
760  mesh,
761  runTime.timeName(),
762  cloud::prefix/cloudDirs[i],
765  false
766  );
767 
768  IOobject* positionsPtr = sprayObjs.lookup
769  (
770  word("positions")
771  );
772 
773  if (positionsPtr)
774  {
775  // Read lagrangian particles
776  // ~~~~~~~~~~~~~~~~~~~~~~~~~
777 
778  Info<< "Identified lagrangian data set: "
779  << cloudDirs[i] << endl;
780 
781  lagrangianPositions.set
782  (
783  cloudI,
784  new Cloud<indexedParticle>
785  (
786  mesh,
787  cloudDirs[i],
788  false
789  )
790  );
791 
792 
793  // Sort particles per cell
794  // ~~~~~~~~~~~~~~~~~~~~~~~
795 
796  cellParticles.set
797  (
798  cloudI,
799  new List<SLList<indexedParticle*>*>
800  (
801  mesh.nCells(),
802  static_cast<SLList<indexedParticle*>*>(nullptr)
803  )
804  );
805 
806  label i = 0;
807 
808  forAllIter
809  (
810  Cloud<indexedParticle>,
811  lagrangianPositions[cloudI],
812  iter
813  )
814  {
815  iter().index() = i++;
816 
817  label celli = iter().cell();
818 
819  // Check
820  if (celli < 0 || celli >= mesh.nCells())
821  {
823  << "Illegal cell number " << celli
824  << " for particle with index "
825  << iter().index()
826  << " at position "
827  << iter().position() << nl
828  << "Cell number should be between 0 and "
829  << mesh.nCells()-1 << nl
830  << "On this mesh the particle should"
831  << " be in cell "
832  << mesh.findCell(iter().position())
833  << exit(FatalError);
834  }
835 
836  if (!cellParticles[cloudI][celli])
837  {
838  cellParticles[cloudI][celli] =
839  new SLList<indexedParticle*>();
840  }
841 
842  cellParticles[cloudI][celli]->append(&iter());
843  }
844 
845  // Read fields
846  // ~~~~~~~~~~~
847 
848  IOobjectList lagrangianObjects
849  (
850  mesh,
851  runTime.timeName(),
852  cloud::prefix/cloudDirs[cloudI],
855  false
856  );
857 
859  (
860  cloudI,
861  lagrangianObjects,
862  lagrangianLabelFields
863  );
864 
866  (
867  cloudI,
868  lagrangianObjects,
869  lagrangianLabelFieldFields
870  );
871 
873  (
874  cloudI,
875  lagrangianObjects,
876  lagrangianScalarFields
877  );
878 
880  (
881  cloudI,
882  lagrangianObjects,
883  lagrangianScalarFieldFields
884  );
885 
887  (
888  cloudI,
889  lagrangianObjects,
890  lagrangianVectorFields
891  );
892 
894  (
895  cloudI,
896  lagrangianObjects,
897  lagrangianVectorFieldFields
898  );
899 
901  (
902  cloudI,
903  lagrangianObjects,
904  lagrangianSphericalTensorFields
905  );
906 
908  (
909  cloudI,
910  lagrangianObjects,
911  lagrangianSphericalTensorFieldFields
912  );
913 
915  (
916  cloudI,
917  lagrangianObjects,
918  lagrangianSymmTensorFields
919  );
920 
922  (
923  cloudI,
924  lagrangianObjects,
925  lagrangianSymmTensorFieldFields
926  );
927 
929  (
930  cloudI,
931  lagrangianObjects,
932  lagrangianTensorFields
933  );
934 
936  (
937  cloudI,
938  lagrangianObjects,
939  lagrangianTensorFieldFields
940  );
941 
942  cloudI++;
943  }
944  }
945 
946  lagrangianPositions.setSize(cloudI);
947  cellParticles.setSize(cloudI);
948  lagrangianLabelFields.setSize(cloudI);
949  lagrangianLabelFieldFields.setSize(cloudI);
950  lagrangianScalarFields.setSize(cloudI);
951  lagrangianScalarFieldFields.setSize(cloudI);
952  lagrangianVectorFields.setSize(cloudI);
953  lagrangianVectorFieldFields.setSize(cloudI);
954  lagrangianSphericalTensorFields.setSize(cloudI);
955  lagrangianSphericalTensorFieldFields.setSize(cloudI);
956  lagrangianSymmTensorFields.setSize(cloudI);
957  lagrangianSymmTensorFieldFields.setSize(cloudI);
958  lagrangianTensorFields.setSize(cloudI);
959  lagrangianTensorFieldFields.setSize(cloudI);
960 
961  Info<< endl;
962 
963  // split the fields over processors
964  for (label proci = 0; proci < mesh.nProcs(); proci++)
965  {
966  Info<< "Processor " << proci << ": field transfer" << endl;
967 
968 
969  // open the database
970  if (!processorDbList.set(proci))
971  {
972  processorDbList.set
973  (
974  proci,
975  new Time
976  (
978  args.rootPath(),
979  args.caseName()
980  /fileName(word("processor") + name(proci))
981  )
982  );
983  }
984  Time& processorDb = processorDbList[proci];
985 
986 
987  processorDb.setTime(runTime);
988 
989  // read the mesh
990  if (!procMeshList.set(proci))
991  {
992  procMeshList.set
993  (
994  proci,
995  new fvMesh
996  (
997  IOobject
998  (
999  regionName,
1000  processorDb.timeName(),
1001  processorDb
1002  )
1003  )
1004  );
1005  }
1006  const fvMesh& procMesh = procMeshList[proci];
1007 
1008  const labelIOList& faceProcAddressing = procAddressing
1009  (
1010  procMeshList,
1011  proci,
1012  "faceProcAddressing",
1013  faceProcAddressingList
1014  );
1015 
1016  const labelIOList& cellProcAddressing = procAddressing
1017  (
1018  procMeshList,
1019  proci,
1020  "cellProcAddressing",
1021  cellProcAddressingList
1022  );
1023 
1024  const labelIOList& boundaryProcAddressing = procAddressing
1025  (
1026  procMeshList,
1027  proci,
1028  "boundaryProcAddressing",
1029  boundaryProcAddressingList
1030  );
1031 
1032 
1033  // FV fields
1034  {
1035  if (!fieldDecomposerList.set(proci))
1036  {
1037  fieldDecomposerList.set
1038  (
1039  proci,
1040  new fvFieldDecomposer
1041  (
1042  mesh,
1043  procMesh,
1044  faceProcAddressing,
1045  cellProcAddressing,
1046  boundaryProcAddressing
1047  )
1048  );
1049  }
1050  const fvFieldDecomposer& fieldDecomposer =
1051  fieldDecomposerList[proci];
1052 
1053  fieldDecomposer.decomposeFields(volScalarFields);
1054  fieldDecomposer.decomposeFields(volVectorFields);
1055  fieldDecomposer.decomposeFields
1056  (
1057  volSphericalTensorFields
1058  );
1059  fieldDecomposer.decomposeFields(volSymmTensorFields);
1060  fieldDecomposer.decomposeFields(volTensorFields);
1061 
1062  fieldDecomposer.decomposeFields(surfaceScalarFields);
1063  fieldDecomposer.decomposeFields(surfaceVectorFields);
1064  fieldDecomposer.decomposeFields
1065  (
1066  surfaceSphericalTensorFields
1067  );
1068  fieldDecomposer.decomposeFields
1069  (
1070  surfaceSymmTensorFields
1071  );
1072  fieldDecomposer.decomposeFields(surfaceTensorFields);
1073 
1074  if (times.size() == 1)
1075  {
1076  // Clear cached decomposer
1077  fieldDecomposerList.set(proci, nullptr);
1078  }
1079  }
1080 
1081  // Dimensioned fields
1082  {
1083  if (!dimFieldDecomposerList.set(proci))
1084  {
1085  dimFieldDecomposerList.set
1086  (
1087  proci,
1088  new dimFieldDecomposer
1089  (
1090  mesh,
1091  procMesh,
1092  faceProcAddressing,
1093  cellProcAddressing
1094  )
1095  );
1096  }
1097  const dimFieldDecomposer& dimDecomposer =
1098  dimFieldDecomposerList[proci];
1099 
1100  dimDecomposer.decomposeFields(dimScalarFields);
1101  dimDecomposer.decomposeFields(dimVectorFields);
1102  dimDecomposer.decomposeFields(dimSphericalTensorFields);
1103  dimDecomposer.decomposeFields(dimSymmTensorFields);
1104  dimDecomposer.decomposeFields(dimTensorFields);
1105 
1106  if (times.size() == 1)
1107  {
1108  dimFieldDecomposerList.set(proci, nullptr);
1109  }
1110  }
1111 
1112 
1113  // Point fields
1114  if
1115  (
1116  pointScalarFields.size()
1117  || pointVectorFields.size()
1118  || pointSphericalTensorFields.size()
1119  || pointSymmTensorFields.size()
1120  || pointTensorFields.size()
1121  )
1122  {
1123  const labelIOList& pointProcAddressing = procAddressing
1124  (
1125  procMeshList,
1126  proci,
1127  "pointProcAddressing",
1128  pointProcAddressingList
1129  );
1130 
1131  const pointMesh& procPMesh = pointMesh::New(procMesh);
1132 
1133  if (!pointFieldDecomposerList.set(proci))
1134  {
1135  pointFieldDecomposerList.set
1136  (
1137  proci,
1138  new pointFieldDecomposer
1139  (
1140  pMesh,
1141  procPMesh,
1142  pointProcAddressing,
1143  boundaryProcAddressing
1144  )
1145  );
1146  }
1147  const pointFieldDecomposer& pointDecomposer =
1148  pointFieldDecomposerList[proci];
1149 
1150  pointDecomposer.decomposeFields(pointScalarFields);
1151  pointDecomposer.decomposeFields(pointVectorFields);
1152  pointDecomposer.decomposeFields
1153  (
1154  pointSphericalTensorFields
1155  );
1156  pointDecomposer.decomposeFields(pointSymmTensorFields);
1157  pointDecomposer.decomposeFields(pointTensorFields);
1158 
1159 
1160  if (times.size() == 1)
1161  {
1162  pointProcAddressingList.set(proci, nullptr);
1163  pointFieldDecomposerList.set(proci, nullptr);
1164  }
1165  }
1166 
1167 
1168  // If there is lagrangian data write it out
1169  forAll(lagrangianPositions, cloudI)
1170  {
1171  if (lagrangianPositions[cloudI].size())
1172  {
1173  lagrangianFieldDecomposer fieldDecomposer
1174  (
1175  mesh,
1176  procMesh,
1177  faceProcAddressing,
1178  cellProcAddressing,
1179  cloudDirs[cloudI],
1180  lagrangianPositions[cloudI],
1181  cellParticles[cloudI]
1182  );
1183 
1184  // Lagrangian fields
1185  {
1186  fieldDecomposer.decomposeFields
1187  (
1188  cloudDirs[cloudI],
1189  lagrangianLabelFields[cloudI]
1190  );
1191  fieldDecomposer.decomposeFieldFields
1192  (
1193  cloudDirs[cloudI],
1194  lagrangianLabelFieldFields[cloudI]
1195  );
1196  fieldDecomposer.decomposeFields
1197  (
1198  cloudDirs[cloudI],
1199  lagrangianScalarFields[cloudI]
1200  );
1201  fieldDecomposer.decomposeFieldFields
1202  (
1203  cloudDirs[cloudI],
1204  lagrangianScalarFieldFields[cloudI]
1205  );
1206  fieldDecomposer.decomposeFields
1207  (
1208  cloudDirs[cloudI],
1209  lagrangianVectorFields[cloudI]
1210  );
1211  fieldDecomposer.decomposeFieldFields
1212  (
1213  cloudDirs[cloudI],
1214  lagrangianVectorFieldFields[cloudI]
1215  );
1216  fieldDecomposer.decomposeFields
1217  (
1218  cloudDirs[cloudI],
1219  lagrangianSphericalTensorFields[cloudI]
1220  );
1221  fieldDecomposer.decomposeFieldFields
1222  (
1223  cloudDirs[cloudI],
1224  lagrangianSphericalTensorFieldFields[cloudI]
1225  );
1226  fieldDecomposer.decomposeFields
1227  (
1228  cloudDirs[cloudI],
1229  lagrangianSymmTensorFields[cloudI]
1230  );
1231  fieldDecomposer.decomposeFieldFields
1232  (
1233  cloudDirs[cloudI],
1234  lagrangianSymmTensorFieldFields[cloudI]
1235  );
1236  fieldDecomposer.decomposeFields
1237  (
1238  cloudDirs[cloudI],
1239  lagrangianTensorFields[cloudI]
1240  );
1241  fieldDecomposer.decomposeFieldFields
1242  (
1243  cloudDirs[cloudI],
1244  lagrangianTensorFieldFields[cloudI]
1245  );
1246  }
1247  }
1248  }
1249 
1250  // Decompose the "uniform" directory in the time region
1251  // directory
1252  decomposeUniform(copyUniform, mesh, processorDb, regionDir);
1253 
1254  // For the first region of a multi-region case additionally
1255  // decompose the "uniform" directory in the time directory
1256  if (regionNames.size() > 1 && regioni == 0)
1257  {
1258  decomposeUniform(copyUniform, mesh, processorDb);
1259  }
1260 
1261  // We have cached all the constant mesh data for the current
1262  // processor. This is only important if running with
1263  // multiple times, otherwise it is just extra storage.
1264  if (times.size() == 1)
1265  {
1266  boundaryProcAddressingList.set(proci, nullptr);
1267  cellProcAddressingList.set(proci, nullptr);
1268  faceProcAddressingList.set(proci, nullptr);
1269  procMeshList.set(proci, nullptr);
1270  processorDbList.set(proci, nullptr);
1271  }
1272  }
1273  }
1274  }
1275  }
1276 
1277  Info<< "\nEnd\n" << endl;
1278 
1279  return 0;
1280 }
1281 
1282 
1283 // ************************************************************************* //
List< instant > instantList
List of instants.
Definition: instantList.H:42
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool rmDir(const fileName &) const =0
Remove a directory 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
bool isFile(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist as a file in the file system?
Definition: POSIX.C:555
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
static const fileName null
An empty fileName.
Definition: fileName.H:97
bool exists(IOobject &io) const
Does ioobject exist. Is either a directory (empty name()) or.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
Definition: POSIX.C:284
static void noParallel()
Remove the parallel options.
Definition: argList.C:174
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
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...
const word & regionDir(const word &regionName)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
stressControl lookup("compactNormalStress") >> compactNormalStress
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:539
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:127
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:283
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.
fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:662
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:198
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:265
objects
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:62
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
virtual void flush() const
Forcibly wait until all output done. Flush any cached data.
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:360
fileName cwd()
Return current working directory path name.
Definition: POSIX.C:241
messageStream Info
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:117
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:158
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
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:404
wordList selectRegionNames(const argList &args, const Time &runTime)
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.