/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation \\/ M anipulation | Copyright (C) 2016-2017 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . Application decomposePar Group grpParallelUtilities Description Automatically decomposes a mesh and fields of a case for parallel execution of OpenFOAM. Usage \b decomposePar [OPTION] Options: - \par -cellDist Write the cell distribution as a labelList, for use with 'manual' decomposition method or as a volScalarField for post-processing. - \par -region \ Decompose named region. Does not check for existence of processor*. - \par -allRegions Decompose all regions in regionProperties. Does not check for existence of processor*. - \par -copyZero Copy \a 0 directory to processor* rather than decompose the fields. - \par -copyUniform Copy any \a uniform directories too. - \par -constant - \par -time xxx:yyy Override controlDict settings and decompose selected times. Does not re-decompose the mesh i.e. does not handle moving mesh or changing mesh cases. - \par -fields Use existing geometry decomposition and convert fields only. - \par -noSets Skip decomposing cellSets, faceSets, pointSets. - \par -force Remove any existing \a processor subdirectories before decomposing the geometry. - \par -ifRequired Only decompose the geometry if the number of domains has changed from a previous decomposition. No \a processor subdirectories will be removed unless the \a -force option is also specified. This option can be used to avoid redundant geometry decomposition (eg, in scripts), but should be used with caution when the underlying (serial) geometry or the decomposition method etc. have been changed between decompositions. \*---------------------------------------------------------------------------*/ #include "OSspecific.H" #include "fvCFD.H" #include "IOobjectList.H" #include "domainDecomposition.H" #include "labelIOField.H" #include "labelFieldIOField.H" #include "scalarIOField.H" #include "scalarFieldIOField.H" #include "vectorIOField.H" #include "vectorFieldIOField.H" #include "sphericalTensorIOField.H" #include "sphericalTensorFieldIOField.H" #include "symmTensorIOField.H" #include "symmTensorFieldIOField.H" #include "tensorIOField.H" #include "tensorFieldIOField.H" #include "pointFields.H" #include "regionProperties.H" #include "readFields.H" #include "dimFieldDecomposer.H" #include "fvFieldDecomposer.H" #include "pointFieldDecomposer.H" #include "lagrangianFieldDecomposer.H" #include "decompositionModel.H" #include "collatedFileOperation.H" #include "faCFD.H" #include "emptyFaPatch.H" #include "faMeshDecomposition.H" #include "faFieldDecomposer.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { const labelIOList& procAddressing ( const PtrList& procMeshList, const label proci, const word& name, PtrList& procAddressingList ) { const fvMesh& procMesh = procMeshList[proci]; if (!procAddressingList.set(proci)) { procAddressingList.set ( proci, new labelIOList ( IOobject ( name, procMesh.facesInstance(), procMesh.meshSubDir, procMesh, IOobject::MUST_READ, IOobject::NO_WRITE, false ) ) ); } return procAddressingList[proci]; } void decomposeUniform ( const bool copyUniform, const domainDecomposition& mesh, const Time& processorDb, const word& regionDir = word::null ) { const Time& runTime = mesh.time(); // Any uniform data to copy/link? const fileName uniformDir(regionDir/"uniform"); if (fileHandler().isDir(runTime.timePath()/uniformDir)) { Info<< "Detected additional non-decomposed files in " << runTime.timePath()/uniformDir << endl; const fileName timePath = fileHandler().filePath(processorDb.timePath()); // If no fields have been decomposed the destination // directory will not have been created so make sure. mkDir(timePath); if (copyUniform || mesh.distributed()) { if (!fileHandler().exists(timePath/uniformDir)) { fileHandler().cp ( runTime.timePath()/uniformDir, timePath/uniformDir ); } } else { // Link with relative paths string parentPath = string("..")/".."; if (regionDir != word::null) { parentPath = parentPath/".."; } fileName currentDir(cwd()); chDir(timePath); if (!fileHandler().exists(uniformDir)) { fileHandler().ln ( parentPath/runTime.timeName()/uniformDir, uniformDir ); } chDir(currentDir); } } } } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { argList::addNote ( "decompose a mesh and fields of a case for parallel execution" ); argList::noParallel(); argList::addOption ( "decomposeParDict", "file", "read decomposePar dictionary from specified location" ); #include "addRegionOption.H" argList::addBoolOption ( "allRegions", "operate on all regions in regionProperties" ); argList::addBoolOption ( "cellDist", "write cell distribution as a labelList - for use with 'manual' " "decomposition method or as a volScalarField for post-processing." ); argList::addBoolOption ( "copyZero", "Copy \a 0 directory to processor* rather than decompose the fields" ); argList::addBoolOption ( "copyUniform", "copy any uniform/ directories too" ); argList::addBoolOption ( "fields", "use existing geometry decomposition and convert fields only" ); argList::addBoolOption ( "noSets", "skip decomposing cellSets, faceSets, pointSets" ); argList::addBoolOption ( "force", "remove existing processor*/ subdirs before decomposing the geometry" ); argList::addBoolOption ( "ifRequired", "only decompose geometry if the number of domains has changed" ); // Include explicit constant options, have zero from time range timeSelector::addOptions(true, false); #include "setRootCase.H" const bool optRegion = args.found("region"); const bool allRegions = args.found("allRegions"); const bool writeCellDist = args.found("cellDist"); const bool copyZero = args.found("copyZero"); const bool copyUniform = args.found("copyUniform"); const bool decomposeSets = !args.found("noSets"); const bool ifRequiredDecomposition = args.found("ifRequired"); bool decomposeFieldsOnly = args.found("fields"); bool forceOverwrite = args.found("force"); // Set time from database #include "createTime.H" // Allow override of time instantList times = timeSelector::selectIfPresent(runTime, args); // Allow override of decomposeParDict location fileName decompDictFile; args.readIfPresent("decomposeParDict", decompDictFile); wordList regionNames; if (allRegions) { Info<< "Decomposing all regions in regionProperties" << nl << endl; regionProperties rp(runTime); wordHashSet names; forAllConstIters(rp, iter) { names.insertMany(iter.object()); } regionNames = names.sortedToc(); } else { regionNames = {fvMesh::defaultRegion}; args.readIfPresent("region", regionNames[0]); } forAll(regionNames, regioni) { const word& regionName = regionNames[regioni]; const word& regionDir = regionName == fvMesh::defaultRegion ? word::null : regionName; Info<< "\n\nDecomposing mesh " << regionName << nl << endl; // Determine the existing processor count directly label nProcs = fileHandler().nProcs(runTime.path(), regionDir); // Get requested numberOfSubdomains directly from the dictionary. // Note: have no mesh yet so cannot use decompositionModel::New const label nDomains = decompositionMethod::nDomains ( IOdictionary ( decompositionModel::selectIO ( IOobject ( "decomposeParDict", runTime.time().system(), regionDir, // use region if non-standard runTime, IOobject::MUST_READ, IOobject::NO_WRITE, false ), decompDictFile ) ) ); if (decomposeFieldsOnly) { // Sanity check on previously decomposed case if (nProcs != nDomains) { FatalErrorInFunction << "Specified -fields, but the case was decomposed with " << nProcs << " domains" << nl << "instead of " << nDomains << " domains as specified in decomposeParDict" << nl << exit(FatalError); } } else if (nProcs) { bool procDirsProblem = true; if (ifRequiredDecomposition && nProcs == nDomains) { // we can reuse the decomposition decomposeFieldsOnly = true; procDirsProblem = false; forceOverwrite = false; Info<< "Using existing processor directories" << nl; } if (allRegions || optRegion) { procDirsProblem = false; forceOverwrite = false; } if (forceOverwrite) { Info<< "Removing " << nProcs << " existing processor directories" << endl; fileHandler().rmDir ( runTime.path()/word("processors"), true // silent (may not have been collated) ); // remove existing processor dirs // reverse order to avoid gaps if someone interrupts the process for (label proci = nProcs-1; proci >= 0; --proci) { fileName procDir ( runTime.path()/(word("processor") + name(proci)) ); fileHandler().rmDir(procDir); } procDirsProblem = false; } if (procDirsProblem) { FatalErrorInFunction << "Case is already decomposed with " << nProcs << " domains, use the -force option or manually" << nl << "remove processor directories before decomposing. e.g.," << nl << " rm -rf " << runTime.path().c_str() << "/processor*" << nl << exit(FatalError); } } Info<< "Create mesh" << endl; domainDecomposition mesh ( IOobject ( regionName, runTime.timeName(), runTime, IOobject::NO_READ, IOobject::NO_WRITE, false ), decompDictFile ); // Decompose the mesh if (!decomposeFieldsOnly) { // Disable buffering when writing mesh since we need to read // it later on when decomposing the fields float bufSz = fileOperations::collatedFileOperation::maxThreadFileBufferSize; fileOperations::collatedFileOperation::maxThreadFileBufferSize = 0; mesh.decomposeMesh(); mesh.writeDecomposition(decomposeSets); if (writeCellDist) { const labelList& procIds = mesh.cellToProc(); // Write the decomposition as labelList for use with 'manual' // decomposition method. labelIOList cellDecomposition ( IOobject ( "cellDecomposition", mesh.facesInstance(), mesh, IOobject::NO_READ, IOobject::NO_WRITE, false ), procIds ); cellDecomposition.write(); Info<< nl << "Wrote decomposition to " << cellDecomposition.objectPath() << " for use in manual decomposition." << endl; // Write as volScalarField for postprocessing. volScalarField cellDist ( IOobject ( "cellDist", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE ), mesh, dimensionedScalar("cellDist", dimless, -1), zeroGradientFvPatchScalarField::typeName ); forAll(procIds, celli) { cellDist[celli] = procIds[celli]; } cellDist.correctBoundaryConditions(); cellDist.write(); Info<< nl << "Wrote decomposition as volScalarField to " << cellDist.name() << " for use in postprocessing." << endl; } fileOperations::collatedFileOperation::maxThreadFileBufferSize = bufSz; } if (copyZero) { // Copy the 0 directory into each of the processor directories fileName prevTimePath; for (label proci = 0; proci < mesh.nProcs(); ++proci) { Time processorDb ( Time::controlDictName, args.rootPath(), args.caseName()/fileName(word("processor") + name(proci)) ); processorDb.setTime(runTime); if (fileHandler().isDir(runTime.timePath())) { // Get corresponding directory name (to handle processors/) const fileName timePath ( fileHandler().objectPath ( IOobject ( "", processorDb.timeName(), processorDb ), word::null ) ); if (timePath != prevTimePath) { Info<< "Processor " << proci << ": copying " << runTime.timePath() << nl << " to " << timePath << endl; fileHandler().cp(runTime.timePath(), timePath); prevTimePath = timePath; } } } } else { // Decompose the field files // Cached processor meshes and maps. These are only preserved if // running with multiple times. PtrList