Go to the documentation of this file.
73 elems[i] = map[elems[i]];
89 const scalar mergeDist
92 if (fullMatch || masterMesh.
nCells() == 0)
112 const string toProcString(
"to" +
name(procI));
120 forAll(masterPatches, patchI)
122 const polyPatch& pp = masterPatches[patchI];
126 isA<processorPolyPatch>(pp)
128 pp.
name().rfind(toProcString)
129 == (pp.
name().size()-toProcString.size())
136 masterFaces.append(meshFaceI++);
140 masterFaces.shrink();
154 forAll(addPatches, patchI)
156 const polyPatch& pp = addPatches[patchI];
158 if (isA<processorPolyPatch>(pp))
160 bool isConnected =
false;
162 for (
label mergedProcI = 0; mergedProcI < procI; mergedProcI++)
164 const string fromProcString
172 if (pp.
name() == fromProcString)
184 addFaces.append(meshFaceI++);
212 const scalar mergeDist,
228 Info<<
"mergeSharedPoints : detected " << pointToMaster.size()
229 <<
" points that are to be merged." <<
endl;
250 forAll(pointProcAddressing, procI)
252 labelList& constructMap = pointProcAddressing[procI];
256 label oldPointI = constructMap[i];
259 label newPointI = map().reversePointMap()[oldPointI];
263 constructMap[i] = -newPointI-2;
265 else if (newPointI >= 0)
267 constructMap[i] = newPointI;
272 <<
"Problem. oldPointI:" << oldPointI
286 const word& regionDir
295 databases[procI].findInstance
302 if (pointsInstance != databases[procI].
timeName())
305 <<
"Your time was specified as " << databases[procI].timeName()
306 <<
" but there is no polyMesh/points in that time." <<
endl
307 <<
"(there is a points file in " << pointsInstance
309 <<
"Please rerun with the correct time specified"
310 <<
" (through the -constant, -time or -latestTime "
311 <<
"(at your option)."
315 Info<<
"Reading points from "
316 << databases[procI].caseName()
317 <<
" for time = " << databases[procI].timeName()
325 databases[procI].findInstance
372 forAll(cellProcAddressing, procI)
374 const labelList& pCells = cellProcAddressing[procI];
378 cellDecomposition.write();
380 Info<<
nl <<
"Wrote decomposition to "
381 << cellDecomposition.objectPath()
382 <<
" for use in manual decomposition." <<
endl;
392 runTime.
setTime(0, oldIndex+1);
407 zeroGradientFvPatchScalarField::typeName
410 forAll(cellDecomposition, cellI)
412 cellDist[cellI] = cellDecomposition[cellI];
417 Info<<
nl <<
"Wrote decomposition as volScalarField to "
418 << cellDist.name() <<
" for use in postprocessing."
427 int main(
int argc,
char *argv[])
431 "reconstruct a mesh using geometric information only"
442 "specify the merge distance relative to the bounding box size "
448 "do (slower) geometric matching on all boundary faces"
453 "write cell distribution as a labelList - for use with 'manual' "
454 "decomposition method or as a volScalarField for post-processing."
461 Info<<
"This is an experimental tool which tries to merge"
462 <<
" individual processor" <<
nl
463 <<
"meshes back into one master mesh. Use it if the original"
464 <<
" master mesh has" <<
nl
465 <<
"been deleted or if the processor meshes have been modified"
466 <<
" (topology change)." <<
nl
467 <<
"This tool will write the resulting mesh to a new time step"
468 <<
" and construct" <<
nl
469 <<
"xxxxProcAddressing files in the processor meshes so"
470 <<
" reconstructPar can be" <<
nl
471 <<
"used to regenerate the fields on the master mesh." <<
nl
473 <<
"Not well tested & use at your own risk!" <<
nl
495 Info<<
"Merge tolerance : " << mergeTol <<
nl
496 <<
"Write tolerance : " << writeTol <<
endl;
501 <<
"Your current settings specify ASCII writing with "
503 <<
"Your merging tolerance (" << mergeTol <<
") is finer than this."
505 <<
"Please change your writeFormat to binary"
506 <<
" or increase the writePrecision" <<
endl
507 <<
"or adjust the merge tolerance (-mergeTol)."
516 Info<<
"Doing geometric matching on all boundary faces." <<
nl <<
endl;
520 Info<<
"Doing geometric matching on correct procBoundaries only."
521 <<
nl <<
"This assumes a correct decomposition." <<
endl;
542 Info<<
"Found " << nProcs <<
" processor directories" <<
nl <<
endl;
550 Info<<
"Reading database "
570 databases[0].times(),
585 databases[procI].setTime(
timeDirs[timeI], timeI);
590 /databases[0].timeName()
594 if (!
isFile(meshPath/
"faces"))
605 const scalar mergeDist = mergeTol*bb.
mag();
607 Info<<
"Overall mesh bounding box : " << bb <<
nl
608 <<
"Relative tolerance : " << mergeTol <<
nl
609 <<
"Absolute matching distance : " << mergeDist <<
nl
620 label masterInternalFaces;
627 Info<<
"Constructing empty mesh to add to." <<
nl <<
endl;
642 for (
label procI = 0; procI < nProcs; procI++)
644 Info<<
"Reading mesh to add from "
645 << databases[procI].caseName()
646 <<
" for time = " << databases[procI].timeName()
663 boundaryProcAddressing[procI] =
692 for (
label mergedI = 0; mergedI < procI; mergedI++)
694 renumber(map().oldCellMap(), cellProcAddressing[mergedI]);
696 renumber(map().oldPointMap(), pointProcAddressing[mergedI]);
701 boundaryProcAddressing[mergedI]
706 renumber(map().addedCellMap(), cellProcAddressing[procI]);
708 renumber(map().addedPointMap(), pointProcAddressing[procI]);
709 renumber(map().addedPatchMap(), boundaryProcAddressing[procI]);
723 Info<<
"\nWriting merged mesh to "
727 if (!masterMesh.
write())
730 <<
"Failed writing polyMesh."
743 Info<<
"Reconstructing the addressing from the processor meshes"
744 <<
" to the newly reconstructed mesh" <<
nl <<
endl;
748 Info<<
"Reading processor " << procI <<
" mesh from "
749 << databases[procI].caseName() <<
endl;
764 Info<<
"Writing pointProcAddressing to "
765 << databases[procI].caseName()
766 /procMesh.facesInstance()
774 "pointProcAddressing",
775 procMesh.facesInstance(),
782 pointProcAddressing[procI]
788 Info<<
"Writing faceProcAddressing to "
789 << databases[procI].caseName()
790 /procMesh.facesInstance()
798 "faceProcAddressing",
799 procMesh.facesInstance(),
811 forAll(faceProcAddr, procFaceI)
813 label masterFaceI = faceProcAddr[procFaceI];
817 !procMesh.isInternalFace(procFaceI)
818 && masterFaceI < masterInternalFaces
824 label procOwn = procMesh.faceOwner()[procFaceI];
825 label masterOwn = masterOwner[masterFaceI];
827 if (cellProcAddressing[procI][procOwn] == masterOwn)
830 faceProcAddr[procFaceI]++;
835 faceProcAddr[procFaceI] =
836 -1 - faceProcAddr[procFaceI];
842 faceProcAddr[procFaceI]++;
846 faceProcAddr.write();
851 Info<<
"Writing cellProcAddressing to "
852 << databases[procI].caseName()
853 /procMesh.facesInstance()
861 "cellProcAddressing",
862 procMesh.facesInstance(),
869 cellProcAddressing[procI]
876 Info<<
"Writing boundaryProcAddressing to "
877 << databases[procI].caseName()
878 /procMesh.facesInstance()
886 "boundaryProcAddressing",
887 procMesh.facesInstance(),
894 boundaryProcAddressing[procI]
vectorField pointField
pointField is a vectorField.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
scalar mag() const
The magnitude of the bounding box span.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from string.
A class for handling file names.
const point & max() const
Maximum describing the bounding box.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
static void addOption(const word &opt, const string ¶m="", const string &usage="")
Add to an option to validOptions with usage information.
static void addNote(const string &)
Add extra notes for the usage information.
static word defaultRegion
Return the default region name.
#define forAll(list, i)
Loop across all elements in list.
autoPtr< mapPolyMesh > mergeSharedPoints(const scalar mergeDist, polyMesh &mesh, labelListList &pointProcAddressing)
IOstream::streamFormat writeFormat() const
Default write format.
static autoPtr< mapAddedPolyMesh > add(fvMesh &mesh0, const fvMesh &mesh1, const faceCoupleInfo &coupleInfo, const bool validBoundary=true)
Inplace add mesh to fvMesh. Maps all stored fields. Returns map.
static const boundBox invertedBox
A very large inverted boundBox: min/max == +/- VGREAT.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Direct mesh changes based on v1.3 polyTopoChange syntax.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
const fileName & facesInstance() const
Return the current instance directory for faces.
virtual bool write() const
Write using setting from DB.
Extract command arguments and options from the supplied argc and argv parameters.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
IOList< label > labelIOList
Label container classes.
const fileName & rootPath() const
Return root path.
Mesh consisting of general polyhedral cells.
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
virtual bool write() const
Write mesh using IO settings from time.
bool set(const label) const
Is element set.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A patch is a list of labels that address the faces in the global face list.
virtual const labelList & faceOwner() const
Return face owner.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
List< cell > cellList
list of cells
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
const point & min() const
Minimum describing the bounding box.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
autoPtr< faceCoupleInfo > determineCoupledFaces(const bool fullMatch, const label procI, const polyMesh &masterMesh, const polyMesh &meshToAdd, const scalar mergeDist)
label nInternalFaces() const
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
errorManip< error > abort(error &err)
boundBox procBounds(const argList &args, const PtrList< Time > &databases, const word ®ionDir)
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
const double e
Elementary charge.
label start() const
Return start label of this patch in the polyMesh face list.
static Map< label > findSharedPoints(const polyMesh &, const scalar mergeTol)
Find topologically and geometrically shared points.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
static instantList timeDirs
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
fileName path() const
Return path.
static unsigned int defaultPrecision()
Return the default precision.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const word & constant() const
Return constant name.
label size() const
Return the number of elements in the PtrList.
static const word null
An empty word.
bool optionFound(const word &opt) const
Return true if the named option is found.
PtrList< labelIOList > & faceProcAddressing
Container for information needed to couple to meshes. When constructed from two meshes and a geometri...
A bounding box defined in terms of the points at its extremities.
instantList select(const instantList &) const
Select a list of Time values that are within the ranges.
A List with indirect addressing.
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
label timeIndex() const
Return current time index.
static void noParallel()
Remove the parallel options.
static const scalar defaultMergeTol
Generic GeometricField class.
Foam::argList args(argc, argv)
static word controlDictName
The default control dictionary name (normally "controlDict")
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const word & name() const
Return name.
int main(int argc, char *argv[])
static void mergePoints(const polyMesh &, const Map< label > &pointToMaster, polyTopoChange &meshMod)
Helper: Merge points.
word name(const complex &)
Return a string representation of a complex.
void writeCellDistance(Time &runTime, const fvMesh &masterMesh, const labelListList &cellProcAddressing)