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