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