/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2016-2021 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 [OPTIONS]
Options:
- \par -allRegions
Decompose all regions in regionProperties. Does not check for
existence of processor*.
- \par -case \
Specify case directory to use (instead of the cwd).
- \par -cellDist
Write the cell distribution as a labelList, for use with 'manual'
decomposition method and as a VTK or volScalarField for visualization.
- \par -constant
Include the 'constant/' dir in the times list.
- \par -copyUniform
Copy any \a uniform directories too.
- \par -copyZero
Copy \a 0 directory to processor* rather than decompose the fields.
- \par -debug-switch \
Specify the value of a registered debug switch. Default is 1
if the value is omitted. (Can be used multiple times)
- \par -decomposeParDict \
Use specified file for decomposePar dictionary.
- \par -dry-run
Test without writing the decomposition. Changes -cellDist to
only write VTK output.
- \par -fields
Use existing geometry decomposition and convert fields only.
- \par fileHandler \
Override the file handler type.
- \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.
- \par -info-switch \
Specify the value of a registered info switch. Default is 1
if the value is omitted. (Can be used multiple times)
- \par -latestTime
Select the latest time.
- \par -lib \
Additional library or library list to load (can be used multiple times).
- \par -noFunctionObjects
Do not execute function objects.
- \par -noSets
Skip decomposing cellSets, faceSets, pointSets.
- \par -noZero
Exclude the \a 0 dir from the times list.
- \par -opt-switch \
Specify the value of a registered optimisation switch (int/bool).
Default is 1 if the value is omitted. (Can be used multiple times)
- \par -region \
Decompose named region. Does not check for existence of processor*.
- \par -time \
Override controlDict settings and decompose selected times. Does not
re-decompose the mesh i.e. does not handle moving mesh or changing
mesh cases. Eg, ':10,20 40:70 1000:', 'none', etc.
- \par -verbose
Additional verbosity.
- \par -doc
Display documentation in browser.
- \par -doc-source
Display source code in browser.
- \par -help
Display short help and exit.
- \par -help-man
Display full help (manpage format) and exit.
- \par -help-notes
Display help notes (description) and exit.
- \par -help-full
Display full help and exit.
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "timeSelector.H"
#include "OSspecific.H"
#include "IOobjectList.H"
#include "decompositionModel.H"
#include "domainDecomposition.H"
#include "domainDecompositionDryRun.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 "fvFieldDecomposer.H"
#include "pointFieldDecomposer.H"
#include "lagrangianFieldDecomposer.H"
#include "emptyFaPatch.H"
#include "faMeshDecomposition.H"
#include "faFieldDecomposer.H"
// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
namespace Foam
{
// Read proc addressing at specific instance.
// Uses polyMesh/fvMesh meshSubDir by default
autoPtr procAddressing
(
const fvMesh& procMesh,
const word& name,
const word& instance,
const word& local = fvMesh::meshSubDir
)
{
return autoPtr::New
(
IOobject
(
name,
instance,
local,
procMesh,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false // do not register
)
);
}
// Read proc addressing at specific instance.
// Uses the finiteArea meshSubDir
autoPtr faProcAddressing
(
const fvMesh& procMesh,
const word& name,
const word& instance,
const word& local = faMesh::meshSubDir
)
{
return procAddressing(procMesh, name, instance, local);
}
// Return cached or read proc addressing from facesInstance
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,
procAddressing(procMesh, name, procMesh.facesInstance())
);
}
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);
}
}
}
} // End namespace Foam
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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",
"Use specified file for decomposePar dictionary"
);
#include "addAllRegionOptions.H"
argList::addDryRunOption
(
"Test without writing the decomposition. "
"Changes -cellDist to only write VTK output."
);
argList::addVerboseOption
(
"Additional verbosity"
);
argList::addOption
(
"domains",
"N",
"Override numberOfSubdomains (-dry-run only)",
true // Advanced option
);
argList::addOption
(
"method",
"name",
"Override decomposition method (-dry-run only)",
true // Advanced option
);
argList::addBoolOption
(
"no-finite-area",
"Suppress finiteArea mesh/field decomposition",
true // Advanced option
);
argList::addBoolOption
(
"no-lagrangian",
"Suppress lagrangian (cloud) decomposition",
true // Advanced option
);
argList::addBoolOption
(
"cellDist",
"Write cell distribution as a labelList - for use with 'manual' "
"decomposition method and as a volScalarField for visualization."
);
argList::addBoolOption
(
"copyZero",
"Copy 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
(
"no-fields",
"Suppress conversion of fields (volume, finite-area, lagrangian)"
);
argList::addBoolOption
(
"no-sets",
"Skip decomposing cellSets, faceSets, pointSets"
);
argList::addOptionCompat("no-sets", {"noSets", 2106});
argList::addBoolOption
(
"force",
"Remove existing processor*/ subdirs before decomposing the geometry"
);
argList::addBoolOption
(
"ifRequired",
"Only decompose geometry if the number of domains has changed"
);
// Allow explicit -constant, have zero from time range
timeSelector::addOptions(true, false); // constant(true), zero(false)
#include "setRootCase.H"
const bool writeCellDist = args.found("cellDist");
// Most of these are ignored for dry-run (not triggered anywhere)
const bool copyZero = args.found("copyZero");
const bool copyUniform = args.found("copyUniform");
const bool decomposeSets = !args.found("no-sets");
const bool decomposeIfRequired = args.found("ifRequired");
const bool doDecompFields = !args.found("no-fields");
const bool doFiniteArea = !args.found("no-finite-area");
const bool doLagrangian = !args.found("no-lagrangian");
bool decomposeFieldsOnly = args.found("fields");
bool forceOverwrite = args.found("force");
// Set time from database
#include "createTime.H"
// Allow override of time (unless dry-run)
instantList times;
if (args.dryRun())
{
Info<< "\ndry-run: ignoring -copy*, -fields, -force, time selection"
<< nl;
}
else
{
if (decomposeFieldsOnly && !doDecompFields)
{
FatalErrorIn(args.executable())
<< "Options -fields and -no-fields are mutually exclusive"
<< " ... giving up" << nl
<< exit(FatalError);
}
if (!doDecompFields)
{
Info<< "Skip decompose of all fields" << nl;
}
if (!doFiniteArea)
{
Info<< "Skip decompose of finiteArea mesh/fields" << nl;
}
if (!doLagrangian)
{
Info<< "Skip decompose of lagrangian positions/fields" << nl;
}
times = timeSelector::selectIfPresent(runTime, args);
}
// Allow override of decomposeParDict location
fileName decompDictFile(args.get("decomposeParDict", ""));
if (!decompDictFile.empty() && !decompDictFile.isAbsolute())
{
decompDictFile = runTime.globalPath()/decompDictFile;
}
// Get region names
#include "getAllRegionOptions.H"
const bool optRegions =
(regionNames.size() != 1 || regionNames[0] != polyMesh::defaultRegion);
if (regionNames.size() == 1 && regionNames[0] != polyMesh::defaultRegion)
{
Info<< "Using region: " << regionNames[0] << nl << endl;
}
forAll(regionNames, regioni)
{
const word& regionName = regionNames[regioni];
const word& regionDir =
(
regionName == polyMesh::defaultRegion ? word::null : regionName
);
if (args.dryRun())
{
Info<< "dry-run: decomposing mesh " << regionName << nl << nl
<< "Create mesh..." << flush;
domainDecompositionDryRun decompTest
(
IOobject
(
regionName,
runTime.timeName(),
runTime,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
),
decompDictFile,
args.getOrDefault