decomposePar.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 -copyUniform \n
47  Copy any \a uniform directories too.
48 
49  - \par -constant
50 
51  - \par -time xxx:yyy \n
52  Override controlDict settings and decompose selected times. Does not
53  re-decompose the mesh i.e. does not handle moving mesh or changing
54  mesh cases.
55 
56  - \par -fields \n
57  Use existing geometry decomposition and convert fields only.
58 
59  - \par -noSets \n
60  Skip decomposing cellSets, faceSets, pointSets.
61 
62  - \par -force \n
63  Remove any existing \a processor subdirectories before decomposing the
64  geometry.
65 
66  - \par -ifRequired \n
67  Only decompose the geometry if the number of domains has changed from a
68  previous decomposition. No \a processor subdirectories will be removed
69  unless the \a -force option is also specified. This option can be used
70  to avoid redundant geometry decomposition (eg, in scripts), but should
71  be used with caution when the underlying (serial) geometry or the
72  decomposition method etc. have been changed between decompositions.
73 
74 \*---------------------------------------------------------------------------*/
75 
76 #include "OSspecific.H"
77 #include "fvCFD.H"
78 #include "IOobjectList.H"
79 #include "domainDecomposition.H"
80 #include "labelIOField.H"
81 #include "labelFieldIOField.H"
82 #include "scalarIOField.H"
83 #include "scalarFieldIOField.H"
84 #include "vectorIOField.H"
85 #include "vectorFieldIOField.H"
86 #include "sphericalTensorIOField.H"
88 #include "symmTensorIOField.H"
89 #include "symmTensorFieldIOField.H"
90 #include "tensorIOField.H"
91 #include "tensorFieldIOField.H"
92 #include "pointFields.H"
93 #include "regionProperties.H"
94 
95 #include "readFields.H"
96 #include "dimFieldDecomposer.H"
97 #include "fvFieldDecomposer.H"
98 #include "pointFieldDecomposer.H"
100 #include "decompositionModel.H"
101 
102 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103 
104 namespace Foam
105 {
106 
107 const labelIOList& procAddressing
108 (
109  const PtrList<fvMesh>& procMeshList,
110  const label proci,
111  const word& name,
112  PtrList<labelIOList>& procAddressingList
113 )
114 {
115  const fvMesh& procMesh = procMeshList[proci];
116 
117  if (!procAddressingList.set(proci))
118  {
119  procAddressingList.set
120  (
121  proci,
122  new labelIOList
123  (
124  IOobject
125  (
126  name,
127  procMesh.facesInstance(),
128  procMesh.meshSubDir,
129  procMesh,
132  false
133  )
134  )
135  );
136  }
137  return procAddressingList[proci];
138 }
139 
140 
141 void decomposeUniform
142 (
143  const bool copyUniform,
144  const domainDecomposition& mesh,
145  const Time& processorDb,
146  const word& regionDir = word::null
147 )
148 {
149  const Time& runTime = mesh.time();
150 
151  // Any uniform data to copy/link?
152  const fileName uniformDir(regionDir/"uniform");
153 
154  if (isDir(runTime.timePath()/uniformDir))
155  {
156  Info<< "Detected additional non-decomposed files in "
157  << runTime.timePath()/uniformDir
158  << endl;
159 
160  const fileName timePath = processorDb.timePath();
161 
162  if (copyUniform || mesh.distributed())
163  {
164  cp
165  (
166  runTime.timePath()/uniformDir,
167  timePath/uniformDir
168  );
169  }
170  else
171  {
172  // link with relative paths
173  string parentPath = string("..")/"..";
174 
175  if (regionDir != word::null)
176  {
177  parentPath = parentPath/"..";
178  }
179 
180  fileName currentDir(cwd());
181  chDir(timePath);
182  ln
183  (
184  parentPath/runTime.timeName()/uniformDir,
185  uniformDir
186  );
187  chDir(currentDir);
188  }
189  }
190 }
191 
192 }
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 int main(int argc, char *argv[])
198 {
200  (
201  "decompose a mesh and fields of a case for parallel execution"
202  );
203 
205  #include "addRegionOption.H"
207  (
208  "allRegions",
209  "operate on all regions in regionProperties"
210  );
212  (
213  "cellDist",
214  "write cell distribution as a labelList - for use with 'manual' "
215  "decomposition method or as a volScalarField for post-processing."
216  );
218  (
219  "copyUniform",
220  "copy any uniform/ directories too"
221  );
223  (
224  "fields",
225  "use existing geometry decomposition and convert fields only"
226  );
228  (
229  "noSets",
230  "skip decomposing cellSets, faceSets, pointSets"
231  );
233  (
234  "force",
235  "remove existing processor*/ subdirs before decomposing the geometry"
236  );
238  (
239  "ifRequired",
240  "only decompose geometry if the number of domains has changed"
241  );
242 
243  // Include explicit constant options, have zero from time range
244  timeSelector::addOptions(true, false);
245 
246  #include "setRootCase.H"
247 
248  bool allRegions = args.optionFound("allRegions");
249  bool writeCellDist = args.optionFound("cellDist");
250  bool copyUniform = args.optionFound("copyUniform");
251  bool decomposeFieldsOnly = args.optionFound("fields");
252  bool decomposeSets = !args.optionFound("noSets");
253  bool forceOverwrite = args.optionFound("force");
254  bool ifRequiredDecomposition = args.optionFound("ifRequired");
255 
256  // Set time from database
257  #include "createTime.H"
258  // Allow override of time
260 
261 
262  wordList regionNames;
263  wordList regionDirs;
264  if (allRegions)
265  {
266  Info<< "Decomposing all regions in regionProperties" << nl << endl;
267  regionProperties rp(runTime);
268  forAllConstIter(HashTable<wordList>, rp, iter)
269  {
270  const wordList& regions = iter();
271  forAll(regions, i)
272  {
273  if (findIndex(regionNames, regions[i]) == -1)
274  {
275  regionNames.append(regions[i]);
276  }
277  }
278  }
279  regionDirs = regionNames;
280  }
281  else
282  {
283  word regionName;
284  if (args.optionReadIfPresent("region", regionName))
285  {
286  regionNames = wordList(1, regionName);
287  regionDirs = regionNames;
288  }
289  else
290  {
291  regionNames = wordList(1, fvMesh::defaultRegion);
292  regionDirs = wordList(1, word::null);
293  }
294  }
295 
296 
297 
298  forAll(regionNames, regioni)
299  {
300  const word& regionName = regionNames[regioni];
301  const word& regionDir = regionDirs[regioni];
302 
303  Info<< "\n\nDecomposing mesh " << regionName << nl << endl;
304 
305 
306  // determine the existing processor count directly
307  label nProcs = 0;
308  while
309  (
310  isDir
311  (
312  runTime.path()
313  / (word("processor") + name(nProcs))
314  / runTime.constant()
315  / regionDir
317  )
318  )
319  {
320  ++nProcs;
321  }
322 
323  // get requested numberOfSubdomains. Note: have no mesh yet so
324  // cannot use decompositionModel::New
325  const label nDomains = readLabel
326  (
327  IOdictionary
328  (
329  IOobject
330  (
331  "decomposeParDict",
332  runTime.time().system(),
333  regionDir, // use region if non-standard
334  runTime,
337  false
338  )
339  ).lookup("numberOfSubdomains")
340  );
341 
342  if (decomposeFieldsOnly)
343  {
344  // Sanity check on previously decomposed case
345  if (nProcs != nDomains)
346  {
348  << "Specified -fields, but the case was decomposed with "
349  << nProcs << " domains"
350  << nl
351  << "instead of " << nDomains
352  << " domains as specified in decomposeParDict"
353  << nl
354  << exit(FatalError);
355  }
356  }
357  else if (nProcs)
358  {
359  bool procDirsProblem = true;
360 
361  if (ifRequiredDecomposition && nProcs == nDomains)
362  {
363  // we can reuse the decomposition
364  decomposeFieldsOnly = true;
365  procDirsProblem = false;
366  forceOverwrite = false;
367 
368  Info<< "Using existing processor directories" << nl;
369  }
370 
371  if (forceOverwrite)
372  {
373  Info<< "Removing " << nProcs
374  << " existing processor directories" << endl;
375 
376  // remove existing processor dirs
377  // reverse order to avoid gaps if someone interrupts the process
378  for (label proci = nProcs-1; proci >= 0; --proci)
379  {
380  fileName procDir
381  (
382  runTime.path()/(word("processor") + name(proci))
383  );
384 
385  rmDir(procDir);
386  }
387 
388  procDirsProblem = false;
389  }
390 
391  if (procDirsProblem)
392  {
394  << "Case is already decomposed with " << nProcs
395  << " domains, use the -force option or manually" << nl
396  << "remove processor directories before decomposing. e.g.,"
397  << nl
398  << " rm -rf " << runTime.path().c_str() << "/processor*"
399  << nl
400  << exit(FatalError);
401  }
402  }
403 
404  Info<< "Create mesh" << endl;
405  domainDecomposition mesh
406  (
407  IOobject
408  (
409  regionName,
410  runTime.timeName(),
411  runTime,
414  false
415  )
416  );
417 
418  // Decompose the mesh
419  if (!decomposeFieldsOnly)
420  {
421  mesh.decomposeMesh();
422 
423  mesh.writeDecomposition(decomposeSets);
424 
425  if (writeCellDist)
426  {
427  const labelList& procIds = mesh.cellToProc();
428 
429  // Write the decomposition as labelList for use with 'manual'
430  // decomposition method.
431  labelIOList cellDecomposition
432  (
433  IOobject
434  (
435  "cellDecomposition",
436  mesh.facesInstance(),
437  mesh,
440  false
441  ),
442  procIds
443  );
444  cellDecomposition.write();
445 
446  Info<< nl << "Wrote decomposition to "
447  << cellDecomposition.objectPath()
448  << " for use in manual decomposition." << endl;
449 
450  // Write as volScalarField for postprocessing.
451  volScalarField cellDist
452  (
453  IOobject
454  (
455  "cellDist",
456  runTime.timeName(),
457  mesh,
460  ),
461  mesh,
462  dimensionedScalar("cellDist", dimless, 0)
463  );
464 
465  forAll(procIds, celli)
466  {
467  cellDist[celli] = procIds[celli];
468  }
469 
470  cellDist.write();
471 
472  Info<< nl << "Wrote decomposition as volScalarField to "
473  << cellDist.name() << " for use in postprocessing."
474  << endl;
475  }
476  }
477 
478 
479 
480  // Caches
481  // ~~~~~~
482  // Cached processor meshes and maps. These are only preserved if running
483  // with multiple times.
484  PtrList<Time> processorDbList(mesh.nProcs());
485  PtrList<fvMesh> procMeshList(mesh.nProcs());
486  PtrList<labelIOList> faceProcAddressingList(mesh.nProcs());
487  PtrList<labelIOList> cellProcAddressingList(mesh.nProcs());
488  PtrList<labelIOList> boundaryProcAddressingList(mesh.nProcs());
489  PtrList<fvFieldDecomposer> fieldDecomposerList(mesh.nProcs());
490  PtrList<dimFieldDecomposer> dimFieldDecomposerList(mesh.nProcs());
491  PtrList<labelIOList> pointProcAddressingList(mesh.nProcs());
492  PtrList<pointFieldDecomposer> pointFieldDecomposerList(mesh.nProcs());
493 
494 
495 
496  // Loop over all times
497  forAll(times, timeI)
498  {
499  runTime.setTime(times[timeI], timeI);
500 
501  Info<< "Time = " << runTime.timeName() << endl;
502 
503  // Search for list of objects for this time
504  IOobjectList objects(mesh, runTime.timeName());
505 
506 
507  // Construct the vol fields
508  // ~~~~~~~~~~~~~~~~~~~~~~~~
509  PtrList<volScalarField> volScalarFields;
510  readFields(mesh, objects, volScalarFields);
511  PtrList<volVectorField> volVectorFields;
512  readFields(mesh, objects, volVectorFields);
513  PtrList<volSphericalTensorField> volSphericalTensorFields;
514  readFields(mesh, objects, volSphericalTensorFields);
515  PtrList<volSymmTensorField> volSymmTensorFields;
516  readFields(mesh, objects, volSymmTensorFields);
517  PtrList<volTensorField> volTensorFields;
518  readFields(mesh, objects, volTensorFields);
519 
520 
521  // Construct the dimensioned fields
522  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
523  PtrList<DimensionedField<scalar, volMesh>> dimScalarFields;
524  readFields(mesh, objects, dimScalarFields);
525  PtrList<DimensionedField<vector, volMesh>> dimVectorFields;
526  readFields(mesh, objects, dimVectorFields);
527  PtrList<DimensionedField<sphericalTensor, volMesh>>
528  dimSphericalTensorFields;
529  readFields(mesh, objects, dimSphericalTensorFields);
530  PtrList<DimensionedField<symmTensor, volMesh>> dimSymmTensorFields;
531  readFields(mesh, objects, dimSymmTensorFields);
532  PtrList<DimensionedField<tensor, volMesh>> dimTensorFields;
533  readFields(mesh, objects, dimTensorFields);
534 
535 
536  // Construct the surface fields
537  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
538  PtrList<surfaceScalarField> surfaceScalarFields;
539  readFields(mesh, objects, surfaceScalarFields);
540  PtrList<surfaceVectorField> surfaceVectorFields;
541  readFields(mesh, objects, surfaceVectorFields);
542  PtrList<surfaceSphericalTensorField> surfaceSphericalTensorFields;
543  readFields(mesh, objects, surfaceSphericalTensorFields);
544  PtrList<surfaceSymmTensorField> surfaceSymmTensorFields;
545  readFields(mesh, objects, surfaceSymmTensorFields);
546  PtrList<surfaceTensorField> surfaceTensorFields;
547  readFields(mesh, objects, surfaceTensorFields);
548 
549 
550  // Construct the point fields
551  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
552  const pointMesh& pMesh = pointMesh::New(mesh);
553 
554  PtrList<pointScalarField> pointScalarFields;
555  readFields(pMesh, objects, pointScalarFields);
556  PtrList<pointVectorField> pointVectorFields;
557  readFields(pMesh, objects, pointVectorFields);
558  PtrList<pointSphericalTensorField> pointSphericalTensorFields;
559  readFields(pMesh, objects, pointSphericalTensorFields);
560  PtrList<pointSymmTensorField> pointSymmTensorFields;
561  readFields(pMesh, objects, pointSymmTensorFields);
562  PtrList<pointTensorField> pointTensorFields;
563  readFields(pMesh, objects, pointTensorFields);
564 
565 
566  // Construct the Lagrangian fields
567  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
568 
569  fileNameList cloudDirs
570  (
571  readDir(runTime.timePath()/cloud::prefix, fileName::DIRECTORY)
572  );
573 
574  // Particles
575  PtrList<Cloud<indexedParticle>> lagrangianPositions
576  (
577  cloudDirs.size()
578  );
579  // Particles per cell
580  PtrList<List<SLList<indexedParticle*>*>> cellParticles
581  (
582  cloudDirs.size()
583  );
584 
585  PtrList<PtrList<labelIOField>> lagrangianLabelFields
586  (
587  cloudDirs.size()
588  );
589  PtrList<PtrList<labelFieldCompactIOField>>
590  lagrangianLabelFieldFields
591  (
592  cloudDirs.size()
593  );
594  PtrList<PtrList<scalarIOField>> lagrangianScalarFields
595  (
596  cloudDirs.size()
597  );
598  PtrList<PtrList<scalarFieldCompactIOField>>
599  lagrangianScalarFieldFields
600  (
601  cloudDirs.size()
602  );
603  PtrList<PtrList<vectorIOField>> lagrangianVectorFields
604  (
605  cloudDirs.size()
606  );
607  PtrList<PtrList<vectorFieldCompactIOField>>
608  lagrangianVectorFieldFields
609  (
610  cloudDirs.size()
611  );
612  PtrList<PtrList<sphericalTensorIOField>>
613  lagrangianSphericalTensorFields
614  (
615  cloudDirs.size()
616  );
617  PtrList<PtrList<sphericalTensorFieldCompactIOField>>
618  lagrangianSphericalTensorFieldFields(cloudDirs.size());
619  PtrList<PtrList<symmTensorIOField>> lagrangianSymmTensorFields
620  (
621  cloudDirs.size()
622  );
623  PtrList<PtrList<symmTensorFieldCompactIOField>>
624  lagrangianSymmTensorFieldFields
625  (
626  cloudDirs.size()
627  );
628  PtrList<PtrList<tensorIOField>> lagrangianTensorFields
629  (
630  cloudDirs.size()
631  );
632  PtrList<PtrList<tensorFieldCompactIOField>>
633  lagrangianTensorFieldFields
634  (
635  cloudDirs.size()
636  );
637 
638  label cloudI = 0;
639 
640  forAll(cloudDirs, i)
641  {
642  IOobjectList sprayObjs
643  (
644  mesh,
645  runTime.timeName(),
646  cloud::prefix/cloudDirs[i],
649  false
650  );
651 
652  IOobject* positionsPtr = sprayObjs.lookup(word("positions"));
653 
654  if (positionsPtr)
655  {
656  // Read lagrangian particles
657  // ~~~~~~~~~~~~~~~~~~~~~~~~~
658 
659  Info<< "Identified lagrangian data set: " << cloudDirs[i]
660  << endl;
661 
662  lagrangianPositions.set
663  (
664  cloudI,
665  new Cloud<indexedParticle>
666  (
667  mesh,
668  cloudDirs[i],
669  false
670  )
671  );
672 
673 
674  // Sort particles per cell
675  // ~~~~~~~~~~~~~~~~~~~~~~~
676 
677  cellParticles.set
678  (
679  cloudI,
680  new List<SLList<indexedParticle*>*>
681  (
682  mesh.nCells(),
683  static_cast<SLList<indexedParticle*>*>(NULL)
684  )
685  );
686 
687  label i = 0;
688 
689  forAllIter
690  (
691  Cloud<indexedParticle>,
692  lagrangianPositions[cloudI],
693  iter
694  )
695  {
696  iter().index() = i++;
697 
698  label celli = iter().cell();
699 
700  // Check
701  if (celli < 0 || celli >= mesh.nCells())
702  {
704  << "Illegal cell number " << celli
705  << " for particle with index " << iter().index()
706  << " at position " << iter().position() << nl
707  << "Cell number should be between 0 and "
708  << mesh.nCells()-1 << nl
709  << "On this mesh the particle should"
710  << " be in cell "
711  << mesh.findCell(iter().position())
712  << exit(FatalError);
713  }
714 
715  if (!cellParticles[cloudI][celli])
716  {
717  cellParticles[cloudI][celli] =
718  new SLList<indexedParticle*>();
719  }
720 
721  cellParticles[cloudI][celli]->append(&iter());
722  }
723 
724  // Read fields
725  // ~~~~~~~~~~~
726 
727  IOobjectList lagrangianObjects
728  (
729  mesh,
730  runTime.timeName(),
731  cloud::prefix/cloudDirs[cloudI],
734  false
735  );
736 
738  (
739  cloudI,
740  lagrangianObjects,
741  lagrangianLabelFields
742  );
743 
745  (
746  cloudI,
747  lagrangianObjects,
748  lagrangianLabelFieldFields
749  );
750 
752  (
753  cloudI,
754  lagrangianObjects,
755  lagrangianScalarFields
756  );
757 
759  (
760  cloudI,
761  lagrangianObjects,
762  lagrangianScalarFieldFields
763  );
764 
766  (
767  cloudI,
768  lagrangianObjects,
769  lagrangianVectorFields
770  );
771 
773  (
774  cloudI,
775  lagrangianObjects,
776  lagrangianVectorFieldFields
777  );
778 
780  (
781  cloudI,
782  lagrangianObjects,
783  lagrangianSphericalTensorFields
784  );
785 
787  (
788  cloudI,
789  lagrangianObjects,
790  lagrangianSphericalTensorFieldFields
791  );
792 
794  (
795  cloudI,
796  lagrangianObjects,
797  lagrangianSymmTensorFields
798  );
799 
801  (
802  cloudI,
803  lagrangianObjects,
804  lagrangianSymmTensorFieldFields
805  );
806 
808  (
809  cloudI,
810  lagrangianObjects,
811  lagrangianTensorFields
812  );
813 
815  (
816  cloudI,
817  lagrangianObjects,
818  lagrangianTensorFieldFields
819  );
820 
821  cloudI++;
822  }
823  }
824 
825  lagrangianPositions.setSize(cloudI);
826  cellParticles.setSize(cloudI);
827  lagrangianLabelFields.setSize(cloudI);
828  lagrangianLabelFieldFields.setSize(cloudI);
829  lagrangianScalarFields.setSize(cloudI);
830  lagrangianScalarFieldFields.setSize(cloudI);
831  lagrangianVectorFields.setSize(cloudI);
832  lagrangianVectorFieldFields.setSize(cloudI);
833  lagrangianSphericalTensorFields.setSize(cloudI);
834  lagrangianSphericalTensorFieldFields.setSize(cloudI);
835  lagrangianSymmTensorFields.setSize(cloudI);
836  lagrangianSymmTensorFieldFields.setSize(cloudI);
837  lagrangianTensorFields.setSize(cloudI);
838  lagrangianTensorFieldFields.setSize(cloudI);
839 
840  Info<< endl;
841 
842  // split the fields over processors
843  for (label proci = 0; proci < mesh.nProcs(); proci++)
844  {
845  Info<< "Processor " << proci << ": field transfer" << endl;
846 
847 
848  // open the database
849  if (!processorDbList.set(proci))
850  {
851  processorDbList.set
852  (
853  proci,
854  new Time
855  (
857  args.rootPath(),
858  args.caseName()
859  /fileName(word("processor") + name(proci))
860  )
861  );
862  }
863  Time& processorDb = processorDbList[proci];
864 
865 
866  processorDb.setTime(runTime);
867 
868  // read the mesh
869  if (!procMeshList.set(proci))
870  {
871  procMeshList.set
872  (
873  proci,
874  new fvMesh
875  (
876  IOobject
877  (
878  regionName,
879  processorDb.timeName(),
880  processorDb
881  )
882  )
883  );
884  }
885  const fvMesh& procMesh = procMeshList[proci];
886 
887  const labelIOList& faceProcAddressing = procAddressing
888  (
889  procMeshList,
890  proci,
891  "faceProcAddressing",
892  faceProcAddressingList
893  );
894 
895  const labelIOList& cellProcAddressing = procAddressing
896  (
897  procMeshList,
898  proci,
899  "cellProcAddressing",
900  cellProcAddressingList
901  );
902 
903  const labelIOList& boundaryProcAddressing = procAddressing
904  (
905  procMeshList,
906  proci,
907  "boundaryProcAddressing",
908  boundaryProcAddressingList
909  );
910 
911 
912  // FV fields
913  {
914  if (!fieldDecomposerList.set(proci))
915  {
916  fieldDecomposerList.set
917  (
918  proci,
919  new fvFieldDecomposer
920  (
921  mesh,
922  procMesh,
923  faceProcAddressing,
924  cellProcAddressing,
925  boundaryProcAddressing
926  )
927  );
928  }
929  const fvFieldDecomposer& fieldDecomposer =
930  fieldDecomposerList[proci];
931 
932  fieldDecomposer.decomposeFields(volScalarFields);
933  fieldDecomposer.decomposeFields(volVectorFields);
934  fieldDecomposer.decomposeFields(volSphericalTensorFields);
935  fieldDecomposer.decomposeFields(volSymmTensorFields);
936  fieldDecomposer.decomposeFields(volTensorFields);
937 
938  fieldDecomposer.decomposeFields(surfaceScalarFields);
939  fieldDecomposer.decomposeFields(surfaceVectorFields);
940  fieldDecomposer.decomposeFields
941  (
942  surfaceSphericalTensorFields
943  );
944  fieldDecomposer.decomposeFields(surfaceSymmTensorFields);
945  fieldDecomposer.decomposeFields(surfaceTensorFields);
946 
947  if (times.size() == 1)
948  {
949  // Clear cached decomposer
950  fieldDecomposerList.set(proci, NULL);
951  }
952  }
953 
954  // Dimensioned fields
955  {
956  if (!dimFieldDecomposerList.set(proci))
957  {
958  dimFieldDecomposerList.set
959  (
960  proci,
961  new dimFieldDecomposer
962  (
963  mesh,
964  procMesh,
965  faceProcAddressing,
966  cellProcAddressing
967  )
968  );
969  }
970  const dimFieldDecomposer& dimDecomposer =
971  dimFieldDecomposerList[proci];
972 
973  dimDecomposer.decomposeFields(dimScalarFields);
974  dimDecomposer.decomposeFields(dimVectorFields);
975  dimDecomposer.decomposeFields(dimSphericalTensorFields);
976  dimDecomposer.decomposeFields(dimSymmTensorFields);
977  dimDecomposer.decomposeFields(dimTensorFields);
978 
979  if (times.size() == 1)
980  {
981  dimFieldDecomposerList.set(proci, NULL);
982  }
983  }
984 
985 
986  // Point fields
987  if
988  (
989  pointScalarFields.size()
990  || pointVectorFields.size()
991  || pointSphericalTensorFields.size()
992  || pointSymmTensorFields.size()
993  || pointTensorFields.size()
994  )
995  {
996  const labelIOList& pointProcAddressing = procAddressing
997  (
998  procMeshList,
999  proci,
1000  "pointProcAddressing",
1001  pointProcAddressingList
1002  );
1003 
1004  const pointMesh& procPMesh = pointMesh::New(procMesh);
1005 
1006  if (!pointFieldDecomposerList.set(proci))
1007  {
1008  pointFieldDecomposerList.set
1009  (
1010  proci,
1011  new pointFieldDecomposer
1012  (
1013  pMesh,
1014  procPMesh,
1015  pointProcAddressing,
1016  boundaryProcAddressing
1017  )
1018  );
1019  }
1020  const pointFieldDecomposer& pointDecomposer =
1021  pointFieldDecomposerList[proci];
1022 
1023  pointDecomposer.decomposeFields(pointScalarFields);
1024  pointDecomposer.decomposeFields(pointVectorFields);
1025  pointDecomposer.decomposeFields(pointSphericalTensorFields);
1026  pointDecomposer.decomposeFields(pointSymmTensorFields);
1027  pointDecomposer.decomposeFields(pointTensorFields);
1028 
1029 
1030  if (times.size() == 1)
1031  {
1032  pointProcAddressingList.set(proci, NULL);
1033  pointFieldDecomposerList.set(proci, NULL);
1034  }
1035  }
1036 
1037 
1038  // If there is lagrangian data write it out
1039  forAll(lagrangianPositions, cloudI)
1040  {
1041  if (lagrangianPositions[cloudI].size())
1042  {
1043  lagrangianFieldDecomposer fieldDecomposer
1044  (
1045  mesh,
1046  procMesh,
1047  faceProcAddressing,
1048  cellProcAddressing,
1049  cloudDirs[cloudI],
1050  lagrangianPositions[cloudI],
1051  cellParticles[cloudI]
1052  );
1053 
1054  // Lagrangian fields
1055  {
1056  fieldDecomposer.decomposeFields
1057  (
1058  cloudDirs[cloudI],
1059  lagrangianLabelFields[cloudI]
1060  );
1061  fieldDecomposer.decomposeFieldFields
1062  (
1063  cloudDirs[cloudI],
1064  lagrangianLabelFieldFields[cloudI]
1065  );
1066  fieldDecomposer.decomposeFields
1067  (
1068  cloudDirs[cloudI],
1069  lagrangianScalarFields[cloudI]
1070  );
1071  fieldDecomposer.decomposeFieldFields
1072  (
1073  cloudDirs[cloudI],
1074  lagrangianScalarFieldFields[cloudI]
1075  );
1076  fieldDecomposer.decomposeFields
1077  (
1078  cloudDirs[cloudI],
1079  lagrangianVectorFields[cloudI]
1080  );
1081  fieldDecomposer.decomposeFieldFields
1082  (
1083  cloudDirs[cloudI],
1084  lagrangianVectorFieldFields[cloudI]
1085  );
1086  fieldDecomposer.decomposeFields
1087  (
1088  cloudDirs[cloudI],
1089  lagrangianSphericalTensorFields[cloudI]
1090  );
1091  fieldDecomposer.decomposeFieldFields
1092  (
1093  cloudDirs[cloudI],
1094  lagrangianSphericalTensorFieldFields[cloudI]
1095  );
1096  fieldDecomposer.decomposeFields
1097  (
1098  cloudDirs[cloudI],
1099  lagrangianSymmTensorFields[cloudI]
1100  );
1101  fieldDecomposer.decomposeFieldFields
1102  (
1103  cloudDirs[cloudI],
1104  lagrangianSymmTensorFieldFields[cloudI]
1105  );
1106  fieldDecomposer.decomposeFields
1107  (
1108  cloudDirs[cloudI],
1109  lagrangianTensorFields[cloudI]
1110  );
1111  fieldDecomposer.decomposeFieldFields
1112  (
1113  cloudDirs[cloudI],
1114  lagrangianTensorFieldFields[cloudI]
1115  );
1116  }
1117  }
1118  }
1119 
1120  // Decompose the "uniform" directory in the time region
1121  // directory
1122  decomposeUniform(copyUniform, mesh, processorDb, regionDir);
1123 
1124  // For the first region of a multi-region case additionally
1125  // decompose the "uniform" directory in the time directory
1126  if (regionNames.size() > 1 && regioni == 0)
1127  {
1128  decomposeUniform(copyUniform, mesh, processorDb);
1129  }
1130 
1131  // We have cached all the constant mesh data for the current
1132  // processor. This is only important if running with multiple
1133  // times, otherwise it is just extra storage.
1134  if (times.size() == 1)
1135  {
1136  boundaryProcAddressingList.set(proci, NULL);
1137  cellProcAddressingList.set(proci, NULL);
1138  faceProcAddressingList.set(proci, NULL);
1139  procMeshList.set(proci, NULL);
1140  processorDbList.set(proci, NULL);
1141  }
1142  }
1143  }
1144  }
1145 
1146  Info<< "\nEnd\n" << endl;
1147 
1148  return 0;
1149 }
1150 
1151 
1152 // ************************************************************************* //
List< instant > instantList
List of instants.
Definition: instantList.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Foam::word regionName
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
Definition: POSIX.C:289
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:620
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:486
const fileName & rootPath() const
Return root path.
Definition: argListI.H:36
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
static const pointMesh & New(const polyMesh &mesh)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:42
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...
dynamicFvMesh & mesh
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
const Time & time() const
Return time.
Definition: IOobject.C:227
regionProperties rp(runTime)
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:284
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:201
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:720
bool rmDir(const fileName &)
Remove a dirctory and its contents.
Definition: POSIX.C:839
static const char nl
Definition: Ostream.H:262
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< word > wordList
A List of words.
Definition: fileName.H:54
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileName cwd()
Return current working directory path name.
Definition: POSIX.C:246
messageStream Info
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:83
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:124
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
fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:527
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.