1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
28 Automatically decomposes a mesh and fields of a case for parallel
29 execution of OpenFOAM.
33 - decomposePar [OPTION]
36 Write the cell distribution as a labelList, for use with 'manual'
37 decomposition method or as a volScalarField for post-processing.
39 \param -region regionName \n
40 Decompose named region. Does not check for existence of processor*.
42 \param -copyUniform \n
43 Copy any \a uniform directories too.
46 Override controlDict settings and use constant directory.
49 Use existing geometry decomposition and convert fields only.
52 Remove any existing \a processor subdirectories before decomposing the
56 Only decompose the geometry if the number of domains has changed from a
57 previous decomposition. No \a processor subdirectories will be removed
58 unless the \a -force option is also specified. This option can be used
59 to avoid redundant geometry decomposition (eg, in scripts), but should
60 be used with caution when the underlying (serial) geometry or the
61 decomposition method etc. have been changed between decompositions.
63 \*---------------------------------------------------------------------------*/
65 #include "OSspecific.H"
67 #include "IOobjectList.H"
68 #include "domainDecomposition.H"
69 #include "labelIOField.H"
70 #include "labelFieldIOField.H"
71 #include "scalarIOField.H"
72 #include "scalarFieldIOField.H"
73 #include "vectorIOField.H"
74 #include "vectorFieldIOField.H"
75 #include "sphericalTensorIOField.H"
76 #include "sphericalTensorFieldIOField.H"
77 #include "symmTensorIOField.H"
78 #include "symmTensorFieldIOField.H"
79 #include "tensorIOField.H"
80 #include "tensorFieldIOField.H"
81 #include "pointFields.H"
83 #include "readFields.H"
84 #include "dimFieldDecomposer.H"
85 #include "fvFieldDecomposer.H"
86 #include "pointFieldDecomposer.H"
87 #include "lagrangianFieldDecomposer.H"
89 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 int main(int argc, char *argv[])
95 "decompose a mesh and fields of a case for parallel execution"
98 argList::noParallel();
99 #include "addRegionOption.H"
100 argList::addBoolOption
103 "write cell distribution as a labelList - for use with 'manual' "
104 "decomposition method or as a volScalarField for post-processing."
106 argList::addBoolOption
109 "copy any uniform/ directories too"
111 argList::addBoolOption
114 "use existing geometry decomposition and convert fields only"
116 argList::addBoolOption
119 "remove existing processor*/ subdirs before decomposing the geometry"
121 argList::addBoolOption
124 "only decompose geometry if the number of domains has changed"
126 argList::addBoolOption
129 "include the 'constant/' dir in the times list"
132 #include "setRootCase.H"
134 word regionName = fvMesh::defaultRegion;
135 word regionDir = word::null;
137 if (args.optionReadIfPresent("region", regionName))
139 regionDir = regionName;
140 Info<< "Decomposing mesh " << regionName << nl << endl;
143 bool writeCellDist = args.optionFound("cellDist");
144 bool copyUniform = args.optionFound("copyUniform");
145 bool decomposeFieldsOnly = args.optionFound("fields");
146 bool forceOverwrite = args.optionFound("force");
147 bool ifRequiredDecomposition = args.optionFound("ifRequired");
149 #include "createTime.H"
151 // Allow -constant to override controlDict settings.
152 if (args.optionFound("constant"))
154 instantList timeDirs = timeSelector::select0(runTime, args);
155 if (runTime.timeName() != runTime.constant())
157 FatalErrorIn(args.executable())
158 << "No '" << runTime.constant() << "' time present." << endl
159 << "Valid times are " << runTime.times()
165 Info<< "Time = " << runTime.timeName() << endl;
167 // determine the existing processor count directly
174 / (word("processor") + name(nProcs))
177 / polyMesh::meshSubDir
184 // get requested numberOfSubdomains
185 const label nDomains = readLabel
192 runTime.time().system(),
193 regionDir, // use region if non-standard
195 IOobject::MUST_READ_IF_MODIFIED,
199 ).lookup("numberOfSubdomains")
202 if (decomposeFieldsOnly)
204 // Sanity check on previously decomposed case
205 if (nProcs != nDomains)
207 FatalErrorIn(args.executable())
208 << "Specified -fields, but the case was decomposed with "
209 << nProcs << " domains"
211 << "instead of " << nDomains
212 << " domains as specified in decomposeParDict"
219 bool procDirsProblem = true;
221 if (ifRequiredDecomposition && nProcs == nDomains)
223 // we can reuse the decomposition
224 decomposeFieldsOnly = true;
225 procDirsProblem = false;
226 forceOverwrite = false;
228 Info<< "Using existing processor directories" << nl;
233 Info<< "Removing " << nProcs
234 << " existing processor directories" << endl;
236 // remove existing processor dirs
237 // reverse order to avoid gaps if someone interrupts the process
238 for (label procI = nProcs-1; procI >= 0; --procI)
242 runTime.path()/(word("processor") + name(procI))
248 procDirsProblem = false;
253 FatalErrorIn(args.executable())
254 << "Case is already decomposed with " << nProcs
255 << " domains, use the -force option or manually" << nl
256 << "remove processor directories before decomposing. e.g.,"
258 << " rm -rf " << runTime.path().c_str() << "/processor*"
264 Info<< "Create mesh" << endl;
265 domainDecomposition mesh
275 // Decompose the mesh
276 if (!decomposeFieldsOnly)
278 mesh.decomposeMesh();
280 mesh.writeDecomposition();
284 const labelList& procIds = mesh.cellToProc();
286 // Write the decomposition as labelList for use with 'manual'
287 // decomposition method.
288 labelIOList cellDecomposition
293 mesh.facesInstance(),
301 cellDecomposition.write();
303 Info<< nl << "Wrote decomposition to "
304 << cellDecomposition.objectPath()
305 << " for use in manual decomposition." << endl;
307 // Write as volScalarField for postprocessing.
308 volScalarField cellDist
319 dimensionedScalar("cellDist", dimless, 0),
320 zeroGradientFvPatchScalarField::typeName
323 forAll(procIds, celli)
325 cellDist[celli] = procIds[celli];
330 Info<< nl << "Wrote decomposition as volScalarField to "
331 << cellDist.name() << " for use in postprocessing."
337 // Search for list of objects for this time
338 IOobjectList objects(mesh, runTime.timeName());
340 // Construct the vol fields
341 // ~~~~~~~~~~~~~~~~~~~~~~~~
342 PtrList<volScalarField> volScalarFields;
343 readFields(mesh, objects, volScalarFields);
345 PtrList<volVectorField> volVectorFields;
346 readFields(mesh, objects, volVectorFields);
348 PtrList<volSphericalTensorField> volSphericalTensorFields;
349 readFields(mesh, objects, volSphericalTensorFields);
351 PtrList<volSymmTensorField> volSymmTensorFields;
352 readFields(mesh, objects, volSymmTensorFields);
354 PtrList<volTensorField> volTensorFields;
355 readFields(mesh, objects, volTensorFields);
358 // Construct the dimensioned fields
359 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
360 PtrList<DimensionedField<scalar, volMesh> > dimScalarFields;
361 readFields(mesh, objects, dimScalarFields);
363 PtrList<DimensionedField<vector, volMesh> > dimVectorFields;
364 readFields(mesh, objects, dimVectorFields);
366 PtrList<DimensionedField<sphericalTensor, volMesh> >
367 dimSphericalTensorFields;
368 readFields(mesh, objects, dimSphericalTensorFields);
370 PtrList<DimensionedField<symmTensor, volMesh> > dimSymmTensorFields;
371 readFields(mesh, objects, dimSymmTensorFields);
373 PtrList<DimensionedField<tensor, volMesh> > dimTensorFields;
374 readFields(mesh, objects, dimTensorFields);
377 // Construct the surface fields
378 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
379 PtrList<surfaceScalarField> surfaceScalarFields;
380 readFields(mesh, objects, surfaceScalarFields);
381 PtrList<surfaceVectorField> surfaceVectorFields;
382 readFields(mesh, objects, surfaceVectorFields);
383 PtrList<surfaceSphericalTensorField> surfaceSphericalTensorFields;
384 readFields(mesh, objects, surfaceSphericalTensorFields);
385 PtrList<surfaceSymmTensorField> surfaceSymmTensorFields;
386 readFields(mesh, objects, surfaceSymmTensorFields);
387 PtrList<surfaceTensorField> surfaceTensorFields;
388 readFields(mesh, objects, surfaceTensorFields);
391 // Construct the point fields
392 // ~~~~~~~~~~~~~~~~~~~~~~~~~~
393 pointMesh pMesh(mesh);
395 PtrList<pointScalarField> pointScalarFields;
396 readFields(pMesh, objects, pointScalarFields);
398 PtrList<pointVectorField> pointVectorFields;
399 readFields(pMesh, objects, pointVectorFields);
401 PtrList<pointSphericalTensorField> pointSphericalTensorFields;
402 readFields(pMesh, objects, pointSphericalTensorFields);
404 PtrList<pointSymmTensorField> pointSymmTensorFields;
405 readFields(pMesh, objects, pointSymmTensorFields);
407 PtrList<pointTensorField> pointTensorFields;
408 readFields(pMesh, objects, pointTensorFields);
411 // Construct the Lagrangian fields
412 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
414 fileNameList cloudDirs
416 readDir(runTime.timePath()/cloud::prefix, fileName::DIRECTORY)
420 PtrList<Cloud<indexedParticle> > lagrangianPositions(cloudDirs.size());
421 // Particles per cell
422 PtrList< List<SLList<indexedParticle*>*> > cellParticles(cloudDirs.size());
424 PtrList<PtrList<labelIOField> > lagrangianLabelFields(cloudDirs.size());
425 PtrList<PtrList<labelFieldCompactIOField> > lagrangianLabelFieldFields
429 PtrList<PtrList<scalarIOField> > lagrangianScalarFields(cloudDirs.size());
430 PtrList<PtrList<scalarFieldCompactIOField> > lagrangianScalarFieldFields
434 PtrList<PtrList<vectorIOField> > lagrangianVectorFields(cloudDirs.size());
435 PtrList<PtrList<vectorFieldCompactIOField> > lagrangianVectorFieldFields
439 PtrList<PtrList<sphericalTensorIOField> > lagrangianSphericalTensorFields
443 PtrList<PtrList<sphericalTensorFieldCompactIOField> >
444 lagrangianSphericalTensorFieldFields(cloudDirs.size());
445 PtrList<PtrList<symmTensorIOField> > lagrangianSymmTensorFields
449 PtrList<PtrList<symmTensorFieldCompactIOField> >
450 lagrangianSymmTensorFieldFields
454 PtrList<PtrList<tensorIOField> > lagrangianTensorFields
458 PtrList<PtrList<tensorFieldCompactIOField> > lagrangianTensorFieldFields
467 IOobjectList sprayObjs
471 cloud::prefix/cloudDirs[i]
474 IOobject* positionsPtr = sprayObjs.lookup("positions");
478 // Read lagrangian particles
479 // ~~~~~~~~~~~~~~~~~~~~~~~~~
481 Info<< "Identified lagrangian data set: " << cloudDirs[i] << endl;
483 lagrangianPositions.set
486 new Cloud<indexedParticle>
495 // Sort particles per cell
496 // ~~~~~~~~~~~~~~~~~~~~~~~
501 new List<SLList<indexedParticle*>*>
504 static_cast<SLList<indexedParticle*>*>(NULL)
512 Cloud<indexedParticle>,
513 lagrangianPositions[cloudI],
517 iter().index() = i++;
519 label celli = iter().cell();
522 if (celli < 0 || celli >= mesh.nCells())
524 FatalErrorIn(args.executable())
525 << "Illegal cell number " << celli
526 << " for particle with index " << iter().index()
527 << " at position " << iter().position() << nl
528 << "Cell number should be between 0 and "
529 << mesh.nCells()-1 << nl
530 << "On this mesh the particle should be in cell "
531 << mesh.findCell(iter().position())
535 if (!cellParticles[cloudI][celli])
537 cellParticles[cloudI][celli] = new SLList<indexedParticle*>
541 cellParticles[cloudI][celli]->append(&iter());
547 IOobjectList lagrangianObjects
551 cloud::prefix/cloudDirs[cloudI]
554 lagrangianFieldDecomposer::readFields
558 lagrangianLabelFields
561 lagrangianFieldDecomposer::readFieldFields
565 lagrangianLabelFieldFields
568 lagrangianFieldDecomposer::readFields
572 lagrangianScalarFields
575 lagrangianFieldDecomposer::readFieldFields
579 lagrangianScalarFieldFields
582 lagrangianFieldDecomposer::readFields
586 lagrangianVectorFields
589 lagrangianFieldDecomposer::readFieldFields
593 lagrangianVectorFieldFields
596 lagrangianFieldDecomposer::readFields
600 lagrangianSphericalTensorFields
603 lagrangianFieldDecomposer::readFieldFields
607 lagrangianSphericalTensorFieldFields
610 lagrangianFieldDecomposer::readFields
614 lagrangianSymmTensorFields
617 lagrangianFieldDecomposer::readFieldFields
621 lagrangianSymmTensorFieldFields
624 lagrangianFieldDecomposer::readFields
628 lagrangianTensorFields
631 lagrangianFieldDecomposer::readFieldFields
635 lagrangianTensorFieldFields
642 lagrangianPositions.setSize(cloudI);
643 cellParticles.setSize(cloudI);
644 lagrangianLabelFields.setSize(cloudI);
645 lagrangianLabelFieldFields.setSize(cloudI);
646 lagrangianScalarFields.setSize(cloudI);
647 lagrangianScalarFieldFields.setSize(cloudI);
648 lagrangianVectorFields.setSize(cloudI);
649 lagrangianVectorFieldFields.setSize(cloudI);
650 lagrangianSphericalTensorFields.setSize(cloudI);
651 lagrangianSphericalTensorFieldFields.setSize(cloudI);
652 lagrangianSymmTensorFields.setSize(cloudI);
653 lagrangianSymmTensorFieldFields.setSize(cloudI);
654 lagrangianTensorFields.setSize(cloudI);
655 lagrangianTensorFieldFields.setSize(cloudI);
658 // Any uniform data to copy/link?
659 fileName uniformDir("uniform");
661 if (isDir(runTime.timePath()/uniformDir))
663 Info<< "Detected additional non-decomposed files in "
664 << runTime.timePath()/uniformDir
674 // split the fields over processors
675 for (label procI = 0; procI < mesh.nProcs(); procI++)
677 Info<< "Processor " << procI << ": field transfer" << endl;
682 Time::controlDictName,
684 args.caseName()/fileName(word("processor") + name(procI))
687 processorDb.setTime(runTime);
689 // remove files remnants that can cause horrible problems
690 // - mut and nut are used to mark the new turbulence models,
691 // their existence prevents old models from being upgraded
693 fileName timeDir(processorDb.path()/processorDb.timeName());
705 processorDb.timeName(),
710 labelIOList faceProcAddressing
714 "faceProcAddressing",
715 procMesh.facesInstance(),
723 labelIOList cellProcAddressing
727 "cellProcAddressing",
728 procMesh.facesInstance(),
736 labelIOList boundaryProcAddressing
740 "boundaryProcAddressing",
741 procMesh.facesInstance(),
751 fvFieldDecomposer fieldDecomposer
757 boundaryProcAddressing
760 fieldDecomposer.decomposeFields(volScalarFields);
761 fieldDecomposer.decomposeFields(volVectorFields);
762 fieldDecomposer.decomposeFields(volSphericalTensorFields);
763 fieldDecomposer.decomposeFields(volSymmTensorFields);
764 fieldDecomposer.decomposeFields(volTensorFields);
766 fieldDecomposer.decomposeFields(surfaceScalarFields);
767 fieldDecomposer.decomposeFields(surfaceVectorFields);
768 fieldDecomposer.decomposeFields(surfaceSphericalTensorFields);
769 fieldDecomposer.decomposeFields(surfaceSymmTensorFields);
770 fieldDecomposer.decomposeFields(surfaceTensorFields);
773 // Dimensioned fields
775 dimFieldDecomposer fieldDecomposer
783 fieldDecomposer.decomposeFields(dimScalarFields);
784 fieldDecomposer.decomposeFields(dimVectorFields);
785 fieldDecomposer.decomposeFields(dimSphericalTensorFields);
786 fieldDecomposer.decomposeFields(dimSymmTensorFields);
787 fieldDecomposer.decomposeFields(dimTensorFields);
794 pointScalarFields.size()
795 || pointVectorFields.size()
796 || pointSphericalTensorFields.size()
797 || pointSymmTensorFields.size()
798 || pointTensorFields.size()
801 labelIOList pointProcAddressing
805 "pointProcAddressing",
806 procMesh.facesInstance(),
814 pointMesh procPMesh(procMesh);
816 pointFieldDecomposer fieldDecomposer
821 boundaryProcAddressing
824 fieldDecomposer.decomposeFields(pointScalarFields);
825 fieldDecomposer.decomposeFields(pointVectorFields);
826 fieldDecomposer.decomposeFields(pointSphericalTensorFields);
827 fieldDecomposer.decomposeFields(pointSymmTensorFields);
828 fieldDecomposer.decomposeFields(pointTensorFields);
832 // If there is lagrangian data write it out
833 forAll(lagrangianPositions, cloudI)
835 if (lagrangianPositions[cloudI].size())
837 lagrangianFieldDecomposer fieldDecomposer
844 lagrangianPositions[cloudI],
845 cellParticles[cloudI]
850 fieldDecomposer.decomposeFields
853 lagrangianLabelFields[cloudI]
855 fieldDecomposer.decomposeFieldFields
858 lagrangianLabelFieldFields[cloudI]
860 fieldDecomposer.decomposeFields
863 lagrangianScalarFields[cloudI]
865 fieldDecomposer.decomposeFieldFields
868 lagrangianScalarFieldFields[cloudI]
870 fieldDecomposer.decomposeFields
873 lagrangianVectorFields[cloudI]
875 fieldDecomposer.decomposeFieldFields
878 lagrangianVectorFieldFields[cloudI]
880 fieldDecomposer.decomposeFields
883 lagrangianSphericalTensorFields[cloudI]
885 fieldDecomposer.decomposeFieldFields
888 lagrangianSphericalTensorFieldFields[cloudI]
890 fieldDecomposer.decomposeFields
893 lagrangianSymmTensorFields[cloudI]
895 fieldDecomposer.decomposeFieldFields
898 lagrangianSymmTensorFieldFields[cloudI]
900 fieldDecomposer.decomposeFields
903 lagrangianTensorFields[cloudI]
905 fieldDecomposer.decomposeFieldFields
908 lagrangianTensorFieldFields[cloudI]
915 // Any non-decomposed data to copy?
916 if (uniformDir.size())
918 const fileName timePath = processorDb.timePath();
920 if (copyUniform || mesh.distributed())
924 runTime.timePath()/uniformDir,
930 // link with relative paths
931 const string parentPath = string("..")/"..";
933 fileName currentDir(cwd());
937 parentPath/runTime.timeName()/uniformDir,
946 Info<< "\nEnd.\n" << endl;
952 // ************************************************************************* //