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