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-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  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 -dict <filename>
70  Specify alternative dictionary for the decomposition.
71 
72 \*---------------------------------------------------------------------------*/
73 
74 #include "processorRunTimes.H"
75 #include "domainDecomposition.H"
76 #include "decompositionMethod.H"
77 #include "argList.H"
78 #include "timeSelector.H"
79 #include "regionProperties.H"
80 
81 #include "labelIOField.H"
82 #include "labelFieldIOField.H"
83 #include "scalarIOField.H"
84 #include "scalarFieldIOField.H"
85 #include "vectorIOField.H"
86 #include "vectorFieldIOField.H"
87 #include "sphericalTensorIOField.H"
89 #include "symmTensorIOField.H"
90 #include "symmTensorFieldIOField.H"
91 #include "tensorIOField.H"
92 #include "tensorFieldIOField.H"
93 
94 #include "readFields.H"
95 #include "dimFieldDecomposer.H"
96 #include "fvFieldDecomposer.H"
97 #include "pointFieldDecomposer.H"
99 
100 using namespace Foam;
101 
102 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103 
104 namespace Foam
105 {
106 
107 void decomposeUniform
108 (
109  const bool copyUniform,
110  const bool distributeUniform,
111  const Time& runTime,
112  const Time& procRunTime,
113  const word& regionDir = word::null
114 )
115 {
116  const fileName uniformDir(regionDir/"uniform");
117 
118  if (fileHandler().isDir(runTime.timePath()/uniformDir))
119  {
120  Info<< "Detected additional non-decomposed files in "
121  << runTime.timePath()/uniformDir
122  << endl;
123 
124  const fileName timePath =
125  fileHandler().filePath(procRunTime.timePath());
126 
127  if (copyUniform || distributeUniform)
128  {
129  if (!fileHandler().exists(timePath/uniformDir))
130  {
131  fileHandler().cp
132  (
133  runTime.timePath()/uniformDir,
134  timePath/uniformDir
135  );
136  }
137  }
138  else
139  {
140  // link with relative paths
141  string parentPath = string("..")/"..";
142 
143  if (regionDir != word::null)
144  {
145  parentPath = parentPath/"..";
146  }
147 
148  fileName currentDir(cwd());
149  chDir(timePath);
150 
151  if (!fileHandler().exists(uniformDir))
152  {
153  fileHandler().ln
154  (
155  parentPath/runTime.timeName()/uniformDir,
156  uniformDir
157  );
158  }
159  chDir(currentDir);
160  }
161  }
162 }
163 
164 
165 void writeDecomposition(const domainDecomposition& meshes)
166 {
167  const labelList& procIds = meshes.cellToProc();
168 
169  // Write the decomposition as labelList for use with 'manual'
170  // decomposition method.
171  labelIOList cellDecomposition
172  (
173  IOobject
174  (
175  "cellDecomposition",
176  meshes.completeMesh().facesInstance(),
177  meshes.completeMesh(),
180  false
181  ),
182  procIds
183  );
184 
185  cellDecomposition.write();
186 
187  Info<< nl << "Wrote decomposition to "
188  << cellDecomposition.relativeObjectPath()
189  << " for use in manual decomposition." << endl;
190 
191  // Write as volScalarField for postprocessing.
192  volScalarField::Internal cellDist
193  (
194  IOobject
195  (
196  "cellDist",
197  meshes.completeMesh().time().timeName(),
198  meshes.completeMesh(),
201  ),
202  meshes.completeMesh(),
203  dimless,
204  scalarField(scalarList(procIds))
205  );
206 
207  cellDist.write();
208 
209  Info<< nl << "Wrote decomposition as volScalarField to "
210  << cellDist.name() << " for use in postprocessing."
211  << endl;
212 }
213 
214 
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 int main(int argc, char *argv[])
221 {
223  (
224  "decompose a mesh and fields of a case for parallel execution"
225  );
226 
228  #include "addDictOption.H"
229  #include "addRegionOption.H"
230  #include "addAllRegionsOption.H"
232  (
233  "cellDist",
234  "write cell distribution as a labelList - for use with 'manual' "
235  "decomposition method or as a volScalarField for post-processing."
236  );
238  (
239  "copyZero",
240  "Copy \a 0 directory to processor* rather than decompose the fields"
241  );
243  (
244  "copyUniform",
245  "copy any uniform/ directories too"
246  );
248  (
249  "fields",
250  "use existing geometry decomposition and convert fields only"
251  );
253  (
254  "noFields",
255  "opposite of -fields; only decompose geometry"
256  );
258  (
259  "noSets",
260  "skip decomposing cellSets, faceSets, pointSets"
261  );
263  (
264  "force",
265  "remove existing processor*/ subdirs before decomposing the geometry"
266  );
267 
268  // Include explicit constant options, have zero from time range
269  timeSelector::addOptions(true, false);
270 
271  #include "setRootCase.H"
272 
273  bool region = args.optionFound("region");
274  bool writeCellDist = args.optionFound("cellDist");
275  bool copyZero = args.optionFound("copyZero");
276  bool copyUniform = args.optionFound("copyUniform");
277  bool decomposeFieldsOnly = args.optionFound("fields");
278  bool decomposeGeomOnly = args.optionFound("noFields");
279  bool decomposeSets = !args.optionFound("noSets");
280  bool forceOverwrite = args.optionFound("force");
281 
282  if (decomposeGeomOnly)
283  {
284  Info<< "Skipping decomposing fields"
285  << nl << endl;
286 
287  if (decomposeFieldsOnly || copyZero)
288  {
290  << "Cannot combine geometry-only decomposition (-noFields)"
291  << " with field decomposition (-fields or -copyZero)"
292  << exit(FatalError);
293  }
294  }
295 
296  // Set time from database
297  Info<< "Create time\n" << endl;
299 
300  // Allow override of time
301  const instantList times = runTimes.selectComplete(args);
302 
303  // Get region names
304  const wordList regionNames =
305  selectRegionNames(args, runTimes.completeTime());
306 
307  // Handle existing decomposition directories
308  {
309  // Determine the processor count from the directories
310  label nProcs = fileHandler().nProcs(runTimes.completeTime().path());
311 
312  if (forceOverwrite)
313  {
314  if (region)
315  {
317  << "Cannot force the decomposition of a single region"
318  << exit(FatalError);
319  }
320 
321  Info<< "Removing " << nProcs
322  << " existing processor directories" << endl;
323 
324  // Remove existing processors directory
325  fileNameList dirs
326  (
328  (
329  runTimes.completeTime().path(),
331  )
332  );
333  forAllReverse(dirs, diri)
334  {
335  const fileName& d = dirs[diri];
336 
337  // Starts with 'processors'
338  if (d.find("processors") == 0)
339  {
340  if (fileHandler().exists(d))
341  {
342  fileHandler().rmDir(d);
343  }
344  }
345 
346  // Starts with 'processor'
347  if (d.find("processor") == 0)
348  {
349  // Check that integer after processor
350  fileName num(d.substr(9));
351  label proci = -1;
352  if (Foam::read(num.c_str(), proci))
353  {
354  if (fileHandler().exists(d))
355  {
356  fileHandler().rmDir(d);
357  }
358  }
359  }
360  }
361  }
362  else if (nProcs && !region && !decomposeFieldsOnly)
363  {
365  << "Case is already decomposed with " << nProcs
366  << " domains, use the -force option or manually" << nl
367  << "remove processor directories before decomposing. e.g.,"
368  << nl
369  << " rm -rf " << runTimes.completeTime().path().c_str()
370  << "/processor*"
371  << nl
372  << exit(FatalError);
373  }
374  }
375 
376  // Get the decomposition dictionary
377  const dictionary decomposeParDict =
378  decompositionMethod::decomposeParDict(runTimes.completeTime());
379 
380  // Decompose all regions
381  forAll(regionNames, regioni)
382  {
383  const word& regionName = regionNames[regioni];
384  const word& regionDir = Foam::regionDir(regionName);
385 
386  Info<< "\n\nDecomposing mesh " << regionName << nl << endl;
387 
388  // Determine the existing processor count directly
389  const label nProcs =
390  fileHandler().nProcs(runTimes.completeTime().path(), regionDir);
391 
392  // Get requested numberOfSubdomains
393  const label nDomains =
394  decomposeParDict.lookup<label>("numberOfSubdomains");
395 
396  // Give file handler a chance to determine the output directory
397  const_cast<fileOperation&>(fileHandler()).setNProcs(nDomains);
398 
399  // Sanity check number of processors in a previously decomposed case
400  if (decomposeFieldsOnly && nProcs != nDomains)
401  {
403  << "Specified -fields, but the case was decomposed with "
404  << nProcs << " domains" << nl << "instead of " << nDomains
405  << " domains as specified in decomposeParDict" << nl
406  << exit(FatalError);
407  }
408 
409  // Get flag to determine whether or not to distribute uniform data
410  const label distributeUniform =
411  decomposeParDict.lookupOrDefault<bool>("distributed", false);
412 
413  // Create meshes
414  Info<< "Create mesh" << endl;
415  domainDecomposition meshes(runTimes, regionName);
416  meshes.readComplete();
417 
418  // Read or generate a decomposition as necessary
419  if (decomposeFieldsOnly)
420  {
421  meshes.readProcs();
422  if (!copyZero)
423  {
424  meshes.readAddressing();
425  meshes.readUpdate();
426  }
427  }
428  else
429  {
430  meshes.decompose();
431  meshes.writeProcs(decomposeSets);
432  if (writeCellDist) writeDecomposition(meshes);
433  fileHandler().flush();
434  }
435 
436  // Field maps. These are preserved if decomposing multiple times.
437  PtrList<fvFieldDecomposer> fieldDecomposerList
438  (
439  meshes.nProcs()
440  );
441  PtrList<dimFieldDecomposer> dimFieldDecomposerList
442  (
443  meshes.nProcs()
444  );
445  PtrList<pointFieldDecomposer> pointFieldDecomposerList
446  (
447  meshes.nProcs()
448  );
449 
450  // Loop over all times
451  forAll(times, timeI)
452  {
453  // Set the time
454  runTimes.setTime(times[timeI], timeI);
455 
456  Info<< "Time = " << runTimes.completeTime().userTimeName() << endl;
457 
458  // Update the meshes, if necessary
460  if (!decomposeFieldsOnly || !copyZero)
461  {
462  state = meshes.readUpdate();
463  }
464 
465  // Write the mesh out, if necessary
466  if (decomposeFieldsOnly)
467  {
468  // Nothing to do
469  }
470  else if (state == fvMesh::POINTS_MOVED)
471  {
472  meshes.writeProcs(false);
473  }
474  else if
475  (
476  state == fvMesh::TOPO_CHANGE
477  || state == fvMesh::TOPO_PATCH_CHANGE
478  )
479  {
480  meshes.writeProcs(decomposeSets);
481  if (writeCellDist) writeDecomposition(meshes);
482  fileHandler().flush();
483  }
484 
485  // Clear the field maps if there has been topology change
486  if (state >= fvMesh::TOPO_CHANGE)
487  {
488  for (label proci = 0; proci < meshes.nProcs(); proci++)
489  {
490  fieldDecomposerList.set(proci, nullptr);
491  dimFieldDecomposerList.set(proci, nullptr);
492  pointFieldDecomposerList.set(proci, nullptr);
493  }
494  }
495 
496  // Decompose the fields at this time
497  if (decomposeGeomOnly)
498  {
499  // Do nothing
500  }
501  else if (copyZero)
502  {
503  // Copy the field files from the <time> directory to the
504  // processor*/<time> directories without altering them
505  const fileName completeTimePath =
506  runTimes.completeTime().timePath();
507 
508  fileName prevProcTimePath;
509  for (label proci = 0; proci < meshes.nProcs(); proci++)
510  {
511  const Time& procRunTime =
512  meshes.procMeshes()[proci].time();
513 
514  if (fileHandler().isDir(completeTimePath))
515  {
516  const fileName procTimePath
517  (
518  fileHandler().objectPath
519  (
520  IOobject
521  (
522  "",
523  procRunTime.timeName(),
524  procRunTime
525  ),
526  word::null
527  )
528  );
529 
530  if (procTimePath != prevProcTimePath)
531  {
532  Info<< "Processor " << proci
533  << ": copying " << completeTimePath << nl
534  << " to " << procTimePath << endl;
535  fileHandler().cp(completeTimePath, procTimePath);
536  prevProcTimePath = procTimePath;
537  }
538  }
539  }
540  }
541  else
542  {
543  // Decompose the fields
544 
545  // Search for list of objects for this time
547  (
548  meshes.completeMesh(),
549  runTimes.completeTime().timeName()
550  );
551 
552  // Construct the vol fields
553  PtrList<volScalarField> volScalarFields;
554  readFields(meshes.completeMesh(), objects, volScalarFields);
555  PtrList<volVectorField> volVectorFields;
556  readFields(meshes.completeMesh(), objects, volVectorFields);
557  PtrList<volSphericalTensorField> volSphericalTensorFields;
558  readFields
559  (
560  meshes.completeMesh(),
561  objects,
562  volSphericalTensorFields
563  );
564  PtrList<volSymmTensorField> volSymmTensorFields;
565  readFields(meshes.completeMesh(), objects, volSymmTensorFields);
566  PtrList<volTensorField> volTensorFields;
567  readFields(meshes.completeMesh(), objects, volTensorFields);
568 
569  // Construct the dimensioned fields
571  readFields(meshes.completeMesh(), objects, dimScalarFields);
573  readFields(meshes.completeMesh(), objects, dimVectorFields);
575  dimSphericalTensorFields;
576  readFields
577  (
578  meshes.completeMesh(),
579  objects,
580  dimSphericalTensorFields
581  );
583  dimSymmTensorFields;
584  readFields(meshes.completeMesh(), objects, dimSymmTensorFields);
586  readFields(meshes.completeMesh(), objects, dimTensorFields);
587 
588  // Construct the surface fields
589  PtrList<surfaceScalarField> surfaceScalarFields;
590  readFields(meshes.completeMesh(), objects, surfaceScalarFields);
591  PtrList<surfaceVectorField> surfaceVectorFields;
592  readFields(meshes.completeMesh(), objects, surfaceVectorFields);
594  surfaceSphericalTensorFields;
595  readFields
596  (
597  meshes.completeMesh(),
598  objects,
599  surfaceSphericalTensorFields
600  );
601  PtrList<surfaceSymmTensorField> surfaceSymmTensorFields;
602  readFields
603  (
604  meshes.completeMesh(),
605  objects,
606  surfaceSymmTensorFields
607  );
608  PtrList<surfaceTensorField> surfaceTensorFields;
609  readFields(meshes.completeMesh(), objects, surfaceTensorFields);
610 
611  // Construct the point fields
612  const pointMesh& pMesh = pointMesh::New(meshes.completeMesh());
613  PtrList<pointScalarField> pointScalarFields;
614  readFields(pMesh, objects, pointScalarFields);
615  PtrList<pointVectorField> pointVectorFields;
616  readFields(pMesh, objects, pointVectorFields);
617  PtrList<pointSphericalTensorField> pointSphericalTensorFields;
618  readFields(pMesh, objects, pointSphericalTensorFields);
619  PtrList<pointSymmTensorField> pointSymmTensorFields;
620  readFields(pMesh, objects, pointSymmTensorFields);
621  PtrList<pointTensorField> pointTensorFields;
622  readFields(pMesh, objects, pointTensorFields);
623 
624  // Construct the Lagrangian fields
625  fileNameList cloudDirs
626  (
628  (
629  runTimes.completeTime().timePath()/cloud::prefix,
631  )
632  );
634  lagrangianPositions(cloudDirs.size());
636  cellParticles(cloudDirs.size());
638  lagrangianLabelFields(cloudDirs.size());
640  lagrangianLabelFieldFields(cloudDirs.size());
642  lagrangianScalarFields(cloudDirs.size());
644  lagrangianScalarFieldFields(cloudDirs.size());
646  lagrangianVectorFields(cloudDirs.size());
648  lagrangianVectorFieldFields(cloudDirs.size());
650  lagrangianSphericalTensorFields(cloudDirs.size());
652  lagrangianSphericalTensorFieldFields(cloudDirs.size());
654  lagrangianSymmTensorFields(cloudDirs.size());
656  lagrangianSymmTensorFieldFields(cloudDirs.size());
658  lagrangianTensorFields(cloudDirs.size());
660  lagrangianTensorFieldFields(cloudDirs.size());
661 
662  label cloudI = 0;
663 
664  forAll(cloudDirs, i)
665  {
666  IOobjectList sprayObjs
667  (
668  meshes.completeMesh(),
669  runTimes.completeTime().timeName(),
670  cloud::prefix/cloudDirs[i],
673  false
674  );
675 
676  IOobject* positionsPtr = sprayObjs.lookup
677  (
678  word("positions")
679  );
680 
681  if (positionsPtr)
682  {
683  // Read lagrangian particles
684  Info<< "Identified lagrangian data set: "
685  << cloudDirs[i] << endl;
686  lagrangianPositions.set
687  (
688  cloudI,
690  (
691  meshes.completeMesh(),
692  cloudDirs[i],
693  false
694  )
695  );
696 
697  // Sort particles per cell
698  cellParticles.set
699  (
700  cloudI,
702  (
703  meshes.completeMesh().nCells(),
704  static_cast<SLList<indexedParticle*>*>(nullptr)
705  )
706  );
707 
708  // Populate the cloud
709  label index = 0;
710  forAllIter
711  (
713  lagrangianPositions[cloudI],
714  iter
715  )
716  {
717  iter().index() = index ++;
718 
719  label celli = iter().cell();
720 
721  // Check
722  if
723  (
724  celli < 0
725  || celli >= meshes.completeMesh().nCells()
726  )
727  {
729  << "Illegal cell number " << celli
730  << " for particle with index "
731  << iter().index()
732  << " at position "
733  << iter().position() << nl
734  << "Cell number should be between 0 and "
735  << meshes.completeMesh().nCells()-1 << nl
736  << "On this mesh the particle should"
737  << " be in cell "
738  << meshes.completeMesh().findCell
739  (iter().position())
740  << exit(FatalError);
741  }
742 
743  if (!cellParticles[cloudI][celli])
744  {
745  cellParticles[cloudI][celli] =
747  }
748 
749  cellParticles[cloudI][celli]->append(&iter());
750  }
751 
752  // Read fields
753  IOobjectList lagrangianObjects
754  (
755  meshes.completeMesh(),
756  runTimes.completeTime().timeName(),
757  cloud::prefix/cloudDirs[cloudI],
760  false
761  );
763  (
764  cloudI,
765  lagrangianObjects,
766  lagrangianLabelFields
767  );
769  (
770  cloudI,
771  lagrangianObjects,
772  lagrangianLabelFieldFields
773  );
775  (
776  cloudI,
777  lagrangianObjects,
778  lagrangianScalarFields
779  );
781  (
782  cloudI,
783  lagrangianObjects,
784  lagrangianScalarFieldFields
785  );
787  (
788  cloudI,
789  lagrangianObjects,
790  lagrangianVectorFields
791  );
793  (
794  cloudI,
795  lagrangianObjects,
796  lagrangianVectorFieldFields
797  );
799  (
800  cloudI,
801  lagrangianObjects,
802  lagrangianSphericalTensorFields
803  );
805  (
806  cloudI,
807  lagrangianObjects,
808  lagrangianSphericalTensorFieldFields
809  );
811  (
812  cloudI,
813  lagrangianObjects,
814  lagrangianSymmTensorFields
815  );
817  (
818  cloudI,
819  lagrangianObjects,
820  lagrangianSymmTensorFieldFields
821  );
823  (
824  cloudI,
825  lagrangianObjects,
826  lagrangianTensorFields
827  );
829  (
830  cloudI,
831  lagrangianObjects,
832  lagrangianTensorFieldFields
833  );
834 
835  cloudI++;
836  }
837  }
838 
839  lagrangianPositions.setSize(cloudI);
840  cellParticles.setSize(cloudI);
841  lagrangianLabelFields.setSize(cloudI);
842  lagrangianLabelFieldFields.setSize(cloudI);
843  lagrangianScalarFields.setSize(cloudI);
844  lagrangianScalarFieldFields.setSize(cloudI);
845  lagrangianVectorFields.setSize(cloudI);
846  lagrangianVectorFieldFields.setSize(cloudI);
847  lagrangianSphericalTensorFields.setSize(cloudI);
848  lagrangianSphericalTensorFieldFields.setSize(cloudI);
849  lagrangianSymmTensorFields.setSize(cloudI);
850  lagrangianSymmTensorFieldFields.setSize(cloudI);
851  lagrangianTensorFields.setSize(cloudI);
852  lagrangianTensorFieldFields.setSize(cloudI);
853 
854  Info<< endl;
855 
856  // split the fields over processors
857  for (label proci = 0; proci < meshes.nProcs(); proci++)
858  {
859  Info<< "Processor " << proci << ": field transfer" << endl;
860 
861  // FV fields
862  {
863  if (!fieldDecomposerList.set(proci))
864  {
865  fieldDecomposerList.set
866  (
867  proci,
869  (
870  meshes.completeMesh(),
871  meshes.procMeshes()[proci],
872  meshes.procFaceAddressing()[proci],
873  meshes.procCellAddressing()[proci],
874  meshes.procFaceAddressingBf()[proci]
875  )
876  );
877  }
878  const fvFieldDecomposer& fieldDecomposer =
879  fieldDecomposerList[proci];
880 
881  fieldDecomposer.decomposeFields(volScalarFields);
882  fieldDecomposer.decomposeFields(volVectorFields);
883  fieldDecomposer.decomposeFields
884  (
885  volSphericalTensorFields
886  );
887  fieldDecomposer.decomposeFields(volSymmTensorFields);
888  fieldDecomposer.decomposeFields(volTensorFields);
889 
890  fieldDecomposer.decomposeFields(surfaceScalarFields);
891  fieldDecomposer.decomposeFields(surfaceVectorFields);
892  fieldDecomposer.decomposeFields
893  (
894  surfaceSphericalTensorFields
895  );
896  fieldDecomposer.decomposeFields
897  (
898  surfaceSymmTensorFields
899  );
900  fieldDecomposer.decomposeFields(surfaceTensorFields);
901 
902  if (times.size() == 1)
903  {
904  // Clear cached decomposer
905  fieldDecomposerList.set(proci, nullptr);
906  }
907  }
908 
909  // Dimensioned fields
910  {
911  if (!dimFieldDecomposerList.set(proci))
912  {
913  dimFieldDecomposerList.set
914  (
915  proci,
917  (
918  meshes.completeMesh(),
919  meshes.procMeshes()[proci],
920  meshes.procFaceAddressing()[proci],
921  meshes.procCellAddressing()[proci]
922  )
923  );
924  }
925  const dimFieldDecomposer& dimDecomposer =
926  dimFieldDecomposerList[proci];
927 
928  dimDecomposer.decomposeFields(dimScalarFields);
929  dimDecomposer.decomposeFields(dimVectorFields);
930  dimDecomposer.decomposeFields(dimSphericalTensorFields);
931  dimDecomposer.decomposeFields(dimSymmTensorFields);
932  dimDecomposer.decomposeFields(dimTensorFields);
933 
934  if (times.size() == 1)
935  {
936  dimFieldDecomposerList.set(proci, nullptr);
937  }
938  }
939 
940  // Point fields
941  if
942  (
943  pointScalarFields.size()
944  || pointVectorFields.size()
945  || pointSphericalTensorFields.size()
946  || pointSymmTensorFields.size()
947  || pointTensorFields.size()
948  )
949  {
950  const pointMesh& procPMesh =
951  pointMesh::New(meshes.procMeshes()[proci]);
952 
953  if (!pointFieldDecomposerList.set(proci))
954  {
955  pointFieldDecomposerList.set
956  (
957  proci,
959  (
960  pMesh,
961  procPMesh,
962  meshes.procPointAddressing()[proci]
963  )
964  );
965  }
966  const pointFieldDecomposer& pointDecomposer =
967  pointFieldDecomposerList[proci];
968 
969  pointDecomposer.decomposeFields(pointScalarFields);
970  pointDecomposer.decomposeFields(pointVectorFields);
971  pointDecomposer.decomposeFields
972  (
973  pointSphericalTensorFields
974  );
975  pointDecomposer.decomposeFields(pointSymmTensorFields);
976  pointDecomposer.decomposeFields(pointTensorFields);
977 
978  if (times.size() == 1)
979  {
980  pointFieldDecomposerList.set(proci, nullptr);
981  }
982  }
983 
984  // If there is lagrangian data write it out
985  forAll(lagrangianPositions, cloudI)
986  {
987  if (lagrangianPositions[cloudI].size())
988  {
989  lagrangianFieldDecomposer fieldDecomposer
990  (
991  meshes.completeMesh(),
992  meshes.procMeshes()[proci],
993  meshes.procFaceAddressing()[proci],
994  meshes.procCellAddressing()[proci],
995  cloudDirs[cloudI],
996  lagrangianPositions[cloudI],
997  cellParticles[cloudI]
998  );
999 
1000  // Lagrangian fields
1001  {
1002  fieldDecomposer.decomposeFields
1003  (
1004  cloudDirs[cloudI],
1005  lagrangianLabelFields[cloudI]
1006  );
1007  fieldDecomposer.decomposeFieldFields
1008  (
1009  cloudDirs[cloudI],
1010  lagrangianLabelFieldFields[cloudI]
1011  );
1012  fieldDecomposer.decomposeFields
1013  (
1014  cloudDirs[cloudI],
1015  lagrangianScalarFields[cloudI]
1016  );
1017  fieldDecomposer.decomposeFieldFields
1018  (
1019  cloudDirs[cloudI],
1020  lagrangianScalarFieldFields[cloudI]
1021  );
1022  fieldDecomposer.decomposeFields
1023  (
1024  cloudDirs[cloudI],
1025  lagrangianVectorFields[cloudI]
1026  );
1027  fieldDecomposer.decomposeFieldFields
1028  (
1029  cloudDirs[cloudI],
1030  lagrangianVectorFieldFields[cloudI]
1031  );
1032  fieldDecomposer.decomposeFields
1033  (
1034  cloudDirs[cloudI],
1035  lagrangianSphericalTensorFields[cloudI]
1036  );
1037  fieldDecomposer.decomposeFieldFields
1038  (
1039  cloudDirs[cloudI],
1040  lagrangianSphericalTensorFieldFields[cloudI]
1041  );
1042  fieldDecomposer.decomposeFields
1043  (
1044  cloudDirs[cloudI],
1045  lagrangianSymmTensorFields[cloudI]
1046  );
1047  fieldDecomposer.decomposeFieldFields
1048  (
1049  cloudDirs[cloudI],
1050  lagrangianSymmTensorFieldFields[cloudI]
1051  );
1052  fieldDecomposer.decomposeFields
1053  (
1054  cloudDirs[cloudI],
1055  lagrangianTensorFields[cloudI]
1056  );
1057  fieldDecomposer.decomposeFieldFields
1058  (
1059  cloudDirs[cloudI],
1060  lagrangianTensorFieldFields[cloudI]
1061  );
1062  }
1063  }
1064  }
1065 
1066  // Decompose the "uniform" directory in the region time
1067  // directory
1068  decomposeUniform
1069  (
1070  copyUniform,
1071  distributeUniform,
1072  runTimes.completeTime(),
1073  meshes.procMeshes()[proci].time(),
1074  regionDir
1075  );
1076 
1077  // For the first region of a multi-region case additionally
1078  // decompose the "uniform" directory in the no-region time
1079  // directory
1080  if (regionNames.size() > 1 && regioni == 0)
1081  {
1082  decomposeUniform
1083  (
1084  copyUniform,
1085  distributeUniform,
1086  runTimes.completeTime(),
1087  meshes.procMeshes()[proci].time()
1088  );
1089  }
1090  }
1091  }
1092  }
1093  }
1094 
1095  Info<< "\nEnd\n" << endl;
1096 
1097  return 0;
1098 }
1099 
1100 
1101 // ************************************************************************* //
void writeProcs(const bool doSets) const
Write the decomposed meshes and associated data.
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.
A class for handling file names.
Definition: fileName.H:79
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:888
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
fileName timePath() const
Return current time path.
Definition: Time.H:287
labelList cellToProc() const
Return the distribution as an FV field for writing.
Template class for non-intrusive linked lists.
Definition: LList.H:51
Dimensioned field decomposer.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const labelListList & procFaceAddressing() const
Access the labels of faces for each processor (see notes above)
label nCells() const
bool exists(IOobject &io) const
Does ioobject exist. Is either a directory (empty name()) or.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
#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:175
const dimensionSet dimless
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
Automatic domain decomposition class for finite-volume meshes.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
static void readFieldFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< CompactIOField< Field< Type >, Type >> > &lagrangianFields)
void decompose()
Decompose the complete mesh to create the processor meshes and.
const word & regionDir(const word &regionName)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual bool ln(const fileName &src, const fileName &dst) const =0
Create a softlink. dst should not exist. Returns true if.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const labelListList & procCellAddressing() const
Access the labels of cells for each processor.
virtual fileName filePath(const bool globalFile, const IOobject &, const word &typeName) const =0
Search for an object. globalFile : also check undecomposed case.
fvMesh::readUpdateState readUpdate()
...
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
static const word null
An empty word.
Definition: word.H:77
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
const PtrList< surfaceLabelField::Boundary > & procFaceAddressingBf() const
Access the labels of faces for each processor (see notes above)
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:207
static dictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
static const char nl
Definition: Ostream.H:260
Point field decomposer.
objects
Finite Volume volume and surface field decomposer.
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:62
void append(const T &a)
Add at tail of list.
Definition: LList.H:172
virtual void flush() const
Forcibly wait until all output done. Flush any cached data.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose llist of fields.
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
fileName cwd()
Return current working directory path name.
Definition: POSIX.C:241
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1647
messageStream Info
const labelListList & procPointAddressing() const
Access the labels of points for each processor.
const PtrList< fvMesh > & procMeshes() const
Access the processor meshes.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
virtual bool write(const bool write=true) const
Write using setting from DB.
static void readFields(const label cloudI, const IOobjectList &lagrangianObjects, PtrList< PtrList< IOField< Type >>> &lagrangianFields)
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose a list of fields.
Foam::argList args(argc, argv)
const fvMesh & completeMesh() const
Access the global mesh.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
A class for handling character strings derived from std::string.
Definition: string.H:76
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
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
Lagrangian field decomposer.
label nProcs() const
Return the number of processors in the decomposition.
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose a list of fields.