102 void decomposeUniform
104 const bool copyUniform,
109 const fileName uniformDir(regionDir/
"uniform");
128 procTimePath/uniformDir
135 string parentPath =
string(
"..")/
"..";
139 parentPath = parentPath/
"..";
181 Info<<
"Wrote decomposition as volScalarField::Internal to "
182 << cellProc.name() <<
" for use in postprocessing"
208 if (!dnl.first_) os <<
nl;
219 int main(
int argc,
char *argv[])
223 "decompose a mesh and fields of a case for parallel execution"
233 "write cell processor indices as a volScalarField::Internal for "
239 "Copy \a 0 directory to processor* rather than decompose the fields"
244 "copy any uniform/ directories too"
249 "use existing geometry decomposition and convert fields only"
254 "opposite of -fields; only decompose geometry"
259 "skip decomposing cellSets, faceSets, pointSets"
264 "remove existing processor*/ subdirs before decomposing the geometry"
282 if (decomposeGeomOnly)
284 Info<<
"Skipping decomposing fields" <<
nl <<
endl;
286 if (decomposeFieldsOnly || copyZero)
289 <<
"Cannot combine geometry-only decomposition (-noFields)"
290 <<
" with field decomposition (-fields or -copyZero)"
311 <<
"Cannot force the decomposition of a single region"
315 const label nProcs0 =
318 Info<<
"Removing " << nProcs0
319 <<
" existing processor directories" <<
nl <<
endl;
335 if (d.find(
"processors") == 0)
344 if (d.find(
"processor") == 0)
366 const label nProcs0 =
369 if (nProcs0 && nProcs0 != runTimes.
nProcs())
372 <<
"Case is already decomposed with " << nProcs0
373 <<
" domains, use the -force option or manually" <<
nl
374 <<
"remove processor directories before decomposing. e.g.,"
391 const word regionDir =
399 const label nDomains =
400 decomposeParDict.
lookup<
label>(
"numberOfSubdomains");
406 if (decomposeFieldsOnly && nProcs != nDomains)
409 <<
"Specified -fields, but the case was decomposed with "
410 << nProcs <<
" domains" <<
nl <<
"instead of " << nDomains
411 <<
" domains as specified in decomposeParDict" <<
nl
420 !(decomposeFieldsOnly && copyZero)
421 && regionMeshes.readDecompose(decomposeSets)
430 writeDecomposition(regionMeshes[regioni]());
438 const bool distributed =
445 runTimes.
setTime(times[timei], timei);
452 !(decomposeFieldsOnly && copyZero)
453 ? regionMeshes.readUpdateDecompose()
458 if (!decomposeFieldsOnly)
460 regionMeshes.writeProcs(decomposeSets);
468 writeDecomposition(regionMeshes[regioni]());
475 if (decomposeGeomOnly)
488 for (
label proci = 0; proci < runTimes.
nProcs(); proci++)
502 if (procTimePath != prevProcTimePath)
504 Info<<
"Processor " << proci
505 <<
": copying " << completeTimePath <<
nl
506 <<
" to " << procTimePath <<
endl;
508 prevProcTimePath = procTimePath;
523 const word regionDir =
526 const delayedNewLine dnl;
531 regionMeshes[regioni];
536 meshes().completeMesh(),
541 Info<< dnl <<
"Decomposing FV fields" <<
endl;
547 meshes().completeMesh(),
548 meshes().procMeshes(),
549 meshes().procFaceAddressing(),
550 meshes().procCellAddressing(),
551 meshes().procFaceAddressingBf()
554 #define DO_FV_VOL_INTERNAL_FIELDS_TYPE(Type, nullArg) \
555 fvDecomposer.decomposeVolInternalFields<Type> \
558 #undef DO_FV_VOL_INTERNAL_FIELDS_TYPE
560 #define DO_FV_VOL_FIELDS_TYPE(Type, nullArg) \
561 fvDecomposer.decomposeVolFields<Type> \
564 #undef DO_FV_VOL_FIELDS_TYPE
566 #define DO_FV_SURFACE_FIELDS_TYPE(Type, nullArg) \
567 fvDecomposer.decomposeFvSurfaceFields<Type> \
570 #undef DO_FV_SURFACE_FIELDS_TYPE
574 Info<< dnl <<
" (no FV fields)" <<
endl;
579 Info<< dnl <<
"Decomposing point fields" <<
endl;
586 meshes().procMeshes(),
587 meshes().procPointAddressing()
590 #define DO_POINT_FIELDS_TYPE(Type, nullArg) \
591 pointDecomposer.decomposeFields<Type> \
594 #undef DO_POINT_FIELDS_TYPE
598 Info<< dnl <<
" (no point fields)" <<
endl;
622 meshes().completeMesh(),
631 if (cloudObjs.lookup(
word(
"positions")))
633 cloudsObjects.
insert(cloudDirs[i], cloudObjs);
638 if (cloudsObjects.
size())
648 string::validate<word>(iter.key());
652 Info<< dnl <<
"Decomposing lagrangian fields for "
666 meshes().completeMesh(),
667 meshes().procMeshes(),
668 meshes().procFaceAddressing(),
669 meshes().procCellAddressing(),
673 lagrangianDecomposer.decomposePositions();
675 #define DO_CLOUD_FIELDS_TYPE(Type, nullArg) \
676 lagrangianDecomposer.decomposeFields<Type> \
678 DO_CLOUD_FIELDS_TYPE(
label, )
680 #undef DO_CLOUD_FIELDS_TYPE
684 Info<< dnl <<
" (no lagrangian fields)"
711 meshes().completeMesh(),
724 || LagrangianFieldDecomposer::decomposes
739 if (LagrangianObjects.
size())
748 const word LagrangianName =
749 string::validate<word>(iter.key());
751 Info<< dnl <<
"Decomposing Lagrangian fields "
752 <<
"for " << LagrangianName <<
endl;
757 meshes().completeMesh(),
758 meshes().procMeshes(),
759 meshes().procFaceAddressing(),
760 meshes().procCellAddressing(),
764 LagrangianDecomposer.decomposePositions();
768 LagrangianFieldDecomposer::decomposes
774 #define DO_LAGRANGIAN_FIELDS_TYPE( \
776 LagrangianDecomposer \
777 .decomposeFields<GeoField<Type>> \
799 #undef DO_LAGRANGIAN_FIELDS_TYPE
808 Info<< dnl <<
" (no lagrangian fields)"
820 if (haveUniform(runTimes))
822 Info<<
"Distributing uniform files" <<
endl;
826 copyUniform || distributed,
839 const word regionDir =
842 if (haveUniform(runTimes, regionDir))
847 regionMeshes[regioni];
849 Info<<
"Distributing uniform files" <<
endl;
853 copyUniform || distributed,
#define forAll(list, i)
Loop across all elements in list.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
An STL-conforming hash table.
label size() const
Return number of elements in table.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
List of IOobjects with searching and retrieving facilities.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Lagrangian field decomposer.
static const word coordinatesName
Name of the coordinates field.
static const word prefix
Instance prefix.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
fileName path() const
Explicitly inherit path from TimePaths to disambiguate from.
word userTimeName() const
Return current user time name with units.
static word controlDictName
The default control dictionary name (normally "controlDict")
fileName timePath() const
Return current time path.
static void addNote(const string &)
Add extra notes for the usage information.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
bool optionFound(const word &opt) const
Return true if the named option is found.
static void noParallel()
Remove the parallel options.
static IOdictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
const word & name() const
Return const reference to name.
Automatic domain decomposition class for finite-volume meshes.
const labelList & cellProc() const
Return the distribution as an FV field for writing.
const fvMesh & completeMesh() const
Access the global mesh.
A class for handling file names.
static const fileName null
An empty fileName.
virtual bool rmDir(const fileName &) const =0
Remove a directory and its contents.
virtual void flush() const
Forcibly wait until all output done. Flush any cached data.
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
virtual fileName filePath(const bool globalFile, const IOobject &) const =0
Search for an object. globalFile : also check undecomposed case.
virtual bool ln(const fileName &src, const fileName &dst) const =0
Create a softlink. dst should not exist. Returns true if.
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
virtual bool isDir(const fileName &, const bool followLink=true) const =0
Does the name exist as a directory in the file system?
virtual bool mkDir(const fileName &, mode_t=0777) const =0
Make directory.
Finite Volume volume and surface field decomposer.
static bool decomposes(const IOobjectList &objects)
Return whether anything in the object list gets decomposed.
const Time & time() const
Return the top-level database.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Lagrangian field decomposer.
static bool decomposes(const IOobjectList &objects)
Return whether anything in the object list gets decomposed.
static const word prefix
The prefix to local: lagrangian.
static bool decomposes(const IOobjectList &objects)
Return whether anything in the object list gets decomposed.
static word defaultRegion
Return the default region name.
instantList selectComplete(const argList &args)
Select the time.
label nProcs() const
Return the number of processors.
const PtrList< Time > & procTimes() const
Access the processor run times.
void setTime(const instant &inst, const label newIndex)
Set the time.
const Time & completeTime() const
Access the complete run time.
A class for handling character strings derived from std::string.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
A class for handling words, derived from string.
static const word null
An empty word.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
int main(int argc, char *argv[])
#define DO_LAGRANGIAN_FIELDS_TYPE(Type, nullArg)
const fileOperation & fileHandler()
Get current file handler.
errorManipArg< error, int > exit(error &err, const int errNo=1)
fileName cwd()
Return current working directory path name.
List< word > wordList
A List of words.
bool read(const char *, int32_t &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimless
const word & regionName(const solver ®ion)
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
List< scalar > scalarList
A List of scalars.
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true)
Read a directory and return the entries as a string list.
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
#define DO_POINT_FIELDS_TYPE(Type, nullArg)
const Foam::wordList regionNames(args.optionFound("allRegions") ? runTime .controlDict().subDict("regionSolvers").toc() :wordList(1, args.optionFound("region") ? args.optionRead< word >("region") :polyMesh::defaultRegion))
Foam::argList args(argc, argv)
const word cloudName(propsDict.lookup("cloudName"))