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