mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: reconstructParMesh - added support for handling multiple time steps
This commit is contained in:
@ -39,7 +39,8 @@ Description
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "argList.H"
|
||||
#include "Time.H"
|
||||
#include "timeSelector.H"
|
||||
|
||||
#include "IOobjectList.H"
|
||||
#include "labelIOList.H"
|
||||
#include "processorPolyPatch.H"
|
||||
@ -278,6 +279,151 @@ autoPtr<mapPolyMesh> mergeSharedPoints
|
||||
}
|
||||
|
||||
|
||||
boundBox procBounds
|
||||
(
|
||||
const argList& args,
|
||||
const PtrList<Time>& databases,
|
||||
const word& regionDir
|
||||
)
|
||||
{
|
||||
boundBox bb = boundBox::invertedBox;
|
||||
|
||||
forAll(databases, procI)
|
||||
{
|
||||
fileName pointsInstance
|
||||
(
|
||||
databases[procI].findInstance
|
||||
(
|
||||
regionDir/polyMesh::meshSubDir,
|
||||
"points"
|
||||
)
|
||||
);
|
||||
|
||||
if (pointsInstance != databases[procI].timeName())
|
||||
{
|
||||
FatalErrorIn(args.executable())
|
||||
<< "Your time was specified as " << databases[procI].timeName()
|
||||
<< " but there is no polyMesh/points in that time." << endl
|
||||
<< "(there is a points file in " << pointsInstance
|
||||
<< ")" << endl
|
||||
<< "Please rerun with the correct time specified"
|
||||
<< " (through the -constant, -time or -latestTime "
|
||||
<< "(at your option)."
|
||||
<< endl << exit(FatalError);
|
||||
}
|
||||
|
||||
Info<< "Reading points from "
|
||||
<< databases[procI].caseName()
|
||||
<< " for time = " << databases[procI].timeName()
|
||||
<< nl << endl;
|
||||
|
||||
pointIOField points
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"points",
|
||||
databases[procI].findInstance
|
||||
(
|
||||
regionDir/polyMesh::meshSubDir,
|
||||
"points"
|
||||
),
|
||||
regionDir/polyMesh::meshSubDir,
|
||||
databases[procI],
|
||||
IOobject::MUST_READ,
|
||||
IOobject::NO_WRITE,
|
||||
false
|
||||
)
|
||||
);
|
||||
|
||||
boundBox domainBb(points, false);
|
||||
|
||||
bb.min() = min(bb.min(), domainBb.min());
|
||||
bb.max() = max(bb.max(), domainBb.max());
|
||||
}
|
||||
|
||||
return bb;
|
||||
}
|
||||
|
||||
|
||||
void writeCellDistance
|
||||
(
|
||||
Time& runTime,
|
||||
const fvMesh& masterMesh,
|
||||
const labelListList& cellProcAddressing
|
||||
|
||||
)
|
||||
{
|
||||
// Write the decomposition as labelList for use with 'manual'
|
||||
// decomposition method.
|
||||
labelIOList cellDecomposition
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"cellDecomposition",
|
||||
masterMesh.facesInstance(),
|
||||
masterMesh,
|
||||
IOobject::NO_READ,
|
||||
IOobject::NO_WRITE,
|
||||
false
|
||||
),
|
||||
masterMesh.nCells()
|
||||
);
|
||||
|
||||
forAll(cellProcAddressing, procI)
|
||||
{
|
||||
const labelList& pCells = cellProcAddressing[procI];
|
||||
UIndirectList<label>(cellDecomposition, pCells) = procI;
|
||||
}
|
||||
|
||||
cellDecomposition.write();
|
||||
|
||||
Info<< nl << "Wrote decomposition to "
|
||||
<< cellDecomposition.objectPath()
|
||||
<< " for use in manual decomposition." << endl;
|
||||
|
||||
|
||||
// Write as volScalarField for postprocessing. Change time to 0
|
||||
// if was 'constant'
|
||||
{
|
||||
const scalar oldTime = runTime.value();
|
||||
const label oldIndex = runTime.timeIndex();
|
||||
if (runTime.timeName() == runTime.constant() && oldIndex == 0)
|
||||
{
|
||||
runTime.setTime(0, oldIndex+1);
|
||||
}
|
||||
|
||||
volScalarField cellDist
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"cellDist",
|
||||
runTime.timeName(),
|
||||
masterMesh,
|
||||
IOobject::NO_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
masterMesh,
|
||||
dimensionedScalar("cellDist", dimless, 0),
|
||||
zeroGradientFvPatchScalarField::typeName
|
||||
);
|
||||
|
||||
forAll(cellDecomposition, cellI)
|
||||
{
|
||||
cellDist[cellI] = cellDecomposition[cellI];
|
||||
}
|
||||
|
||||
cellDist.write();
|
||||
|
||||
Info<< nl << "Wrote decomposition as volScalarField to "
|
||||
<< cellDist.name() << " for use in postprocessing."
|
||||
<< endl;
|
||||
|
||||
// Restore time
|
||||
runTime.setTime(oldTime, oldIndex);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
int main(int argc, char *argv[])
|
||||
{
|
||||
argList::addNote
|
||||
@ -390,7 +536,7 @@ int main(int argc, char *argv[])
|
||||
Info<< "Found " << nProcs << " processor directories" << nl << endl;
|
||||
|
||||
|
||||
// Read all databases.
|
||||
// Read all time databases
|
||||
PtrList<Time> databases(nProcs);
|
||||
|
||||
forAll(databases, procI)
|
||||
@ -409,90 +555,46 @@ int main(int argc, char *argv[])
|
||||
args.caseName()/fileName(word("processor") + name(procI))
|
||||
)
|
||||
);
|
||||
}
|
||||
|
||||
Time& procTime = databases[procI];
|
||||
|
||||
instantList Times = procTime.times();
|
||||
// use the times list from the master processor
|
||||
// and select a subset based on the command-line options
|
||||
instantList Times = timeSelector::select
|
||||
(
|
||||
databases[0].times(),
|
||||
args
|
||||
);
|
||||
|
||||
// set startTime and endTime depending on -time and -latestTime options
|
||||
#include "checkTimeOptions.H"
|
||||
|
||||
procTime.setTime(Times[startTime], startTime);
|
||||
|
||||
if (procI > 0 && databases[procI-1].value() != procTime.value())
|
||||
if (Times.empty())
|
||||
{
|
||||
FatalErrorIn(args.executable())
|
||||
<< "Time not equal on processors." << nl
|
||||
<< "Processor:" << procI-1
|
||||
<< " time:" << databases[procI-1].value() << nl
|
||||
<< "Processor:" << procI
|
||||
<< " time:" << procTime.value()
|
||||
<< "No times selected"
|
||||
<< exit(FatalError);
|
||||
}
|
||||
|
||||
|
||||
// Loop over all times
|
||||
forAll(Times, timeI)
|
||||
{
|
||||
// Set time for global database
|
||||
runTime.setTime(Times[timeI], timeI);
|
||||
|
||||
Info<< "Time = " << runTime.timeName() << endl << endl;
|
||||
|
||||
// Set time for all databases
|
||||
forAll(databases, procI)
|
||||
{
|
||||
databases[procI].setTime(Times[timeI], timeI);
|
||||
}
|
||||
|
||||
// Set master time
|
||||
Info<< "Setting master time to " << databases[0].timeName() << nl << endl;
|
||||
runTime.setTime(databases[0]);
|
||||
|
||||
|
||||
// Read point on individual processors to determine merge tolerance
|
||||
// (otherwise single cell domains might give problems)
|
||||
|
||||
boundBox bb = boundBox::invertedBox;
|
||||
|
||||
for (label procI = 0; procI < nProcs; procI++)
|
||||
{
|
||||
fileName pointsInstance
|
||||
(
|
||||
databases[procI].findInstance
|
||||
(
|
||||
regionDir/polyMesh::meshSubDir,
|
||||
"points"
|
||||
)
|
||||
);
|
||||
|
||||
if (pointsInstance != databases[procI].timeName())
|
||||
{
|
||||
FatalErrorIn(args.executable())
|
||||
<< "Your time was specified as " << databases[procI].timeName()
|
||||
<< " but there is no polyMesh/points in that time." << endl
|
||||
<< "(there is a points file in " << pointsInstance
|
||||
<< ")" << endl
|
||||
<< "Please rerun with the correct time specified"
|
||||
<< " (through the -constant, -time or -latestTime "
|
||||
<< "(at your option)."
|
||||
<< endl << exit(FatalError);
|
||||
}
|
||||
|
||||
Info<< "Reading points from "
|
||||
<< databases[procI].caseName()
|
||||
<< " for time = " << databases[procI].timeName()
|
||||
<< nl << endl;
|
||||
|
||||
pointIOField points
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"points",
|
||||
databases[procI].findInstance
|
||||
(
|
||||
regionDir/polyMesh::meshSubDir,
|
||||
"points"
|
||||
),
|
||||
regionDir/polyMesh::meshSubDir,
|
||||
databases[procI],
|
||||
IOobject::MUST_READ,
|
||||
IOobject::NO_WRITE,
|
||||
false
|
||||
)
|
||||
);
|
||||
|
||||
boundBox domainBb(points, false);
|
||||
|
||||
bb.min() = min(bb.min(), domainBb.min());
|
||||
bb.max() = max(bb.max(), domainBb.max());
|
||||
}
|
||||
const boundBox bb = procBounds(args, databases, regionDir);
|
||||
const scalar mergeDist = mergeTol*bb.mag();
|
||||
|
||||
Info<< "Overall mesh bounding box : " << bb << nl
|
||||
@ -575,8 +677,8 @@ int main(int argc, char *argv[])
|
||||
couples
|
||||
);
|
||||
|
||||
// Update all addressing so xxProcAddressing points to correct item
|
||||
// in masterMesh.
|
||||
// Update all addressing so xxProcAddressing points to correct
|
||||
// item in masterMesh.
|
||||
|
||||
// Processors that were already in masterMesh
|
||||
for (label mergedI = 0; mergedI < procI; mergedI++)
|
||||
@ -585,7 +687,11 @@ int main(int argc, char *argv[])
|
||||
renumber(map().oldFaceMap(), faceProcAddressing[mergedI]);
|
||||
renumber(map().oldPointMap(), pointProcAddressing[mergedI]);
|
||||
// Note: boundary is special since can contain -1.
|
||||
renumber(map().oldPatchMap(), boundaryProcAddressing[mergedI]);
|
||||
renumber
|
||||
(
|
||||
map().oldPatchMap(),
|
||||
boundaryProcAddressing[mergedI]
|
||||
);
|
||||
}
|
||||
|
||||
// Added processor
|
||||
@ -617,78 +723,9 @@ int main(int argc, char *argv[])
|
||||
<< exit(FatalError);
|
||||
}
|
||||
|
||||
|
||||
|
||||
if (writeCellDist)
|
||||
{
|
||||
// Write the decomposition as labelList for use with 'manual'
|
||||
// decomposition method.
|
||||
labelIOList cellDecomposition
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"cellDecomposition",
|
||||
masterMesh.facesInstance(),
|
||||
masterMesh,
|
||||
IOobject::NO_READ,
|
||||
IOobject::NO_WRITE,
|
||||
false
|
||||
),
|
||||
masterMesh.nCells()
|
||||
);
|
||||
|
||||
forAll(cellProcAddressing, procI)
|
||||
{
|
||||
const labelList& pCells = cellProcAddressing[procI];
|
||||
UIndirectList<label>(cellDecomposition, pCells) = procI;
|
||||
}
|
||||
|
||||
cellDecomposition.write();
|
||||
|
||||
Info<< nl << "Wrote decomposition to "
|
||||
<< cellDecomposition.objectPath()
|
||||
<< " for use in manual decomposition." << endl;
|
||||
|
||||
|
||||
// Write as volScalarField for postprocessing. Change time to 0
|
||||
// if was 'constant'
|
||||
{
|
||||
const scalar oldTime = runTime.value();
|
||||
const label oldIndex = runTime.timeIndex();
|
||||
if (runTime.timeName() == runTime.constant() && oldIndex == 0)
|
||||
{
|
||||
runTime.setTime(0, oldIndex+1);
|
||||
}
|
||||
|
||||
volScalarField cellDist
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"cellDist",
|
||||
runTime.timeName(),
|
||||
masterMesh,
|
||||
IOobject::NO_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
masterMesh,
|
||||
dimensionedScalar("cellDist", dimless, 0),
|
||||
zeroGradientFvPatchScalarField::typeName
|
||||
);
|
||||
|
||||
forAll(cellDecomposition, cellI)
|
||||
{
|
||||
cellDist[cellI] = cellDecomposition[cellI];
|
||||
}
|
||||
|
||||
cellDist.write();
|
||||
|
||||
Info<< nl << "Wrote decomposition as volScalarField to "
|
||||
<< cellDist.name() << " for use in postprocessing."
|
||||
<< endl;
|
||||
|
||||
// Restore time
|
||||
runTime.setTime(oldTime, oldIndex);
|
||||
}
|
||||
writeCellDistance(runTime, masterMesh, cellProcAddressing);
|
||||
}
|
||||
}
|
||||
|
||||
@ -851,6 +888,7 @@ int main(int argc, char *argv[])
|
||||
|
||||
Info<< endl;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
Info<< "End.\n" << endl;
|
||||
|
||||
Reference in New Issue
Block a user