mirror of
https://github.com/OpenFOAM/OpenFOAM-6.git
synced 2025-12-08 06:57:46 +00:00
reconstructParMesh: Added -allRegions option
This commit is contained in:
@ -1,8 +1,10 @@
|
|||||||
EXE_INC = \
|
EXE_INC = \
|
||||||
-I$(LIB_SRC)/dynamicMesh/lnInclude \
|
-I$(LIB_SRC)/dynamicMesh/lnInclude \
|
||||||
-I$(LIB_SRC)/finiteVolume/lnInclude \
|
-I$(LIB_SRC)/finiteVolume/lnInclude \
|
||||||
-I$(LIB_SRC)/meshTools/lnInclude
|
-I$(LIB_SRC)/meshTools/lnInclude \
|
||||||
|
-I$(LIB_SRC)/regionModels/regionModel/lnInclude
|
||||||
|
|
||||||
EXE_LIBS = \
|
EXE_LIBS = \
|
||||||
-ldynamicMesh \
|
-ldynamicMesh \
|
||||||
-lmeshTools
|
-lmeshTools \
|
||||||
|
-lregionModels
|
||||||
|
|||||||
@ -2,7 +2,7 @@
|
|||||||
========= |
|
========= |
|
||||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
\\ / O peration |
|
\\ / O peration |
|
||||||
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
|
\\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
|
||||||
\\/ M anipulation |
|
\\/ M anipulation |
|
||||||
-------------------------------------------------------------------------------
|
-------------------------------------------------------------------------------
|
||||||
License
|
License
|
||||||
@ -50,6 +50,7 @@ Description
|
|||||||
#include "fvMeshAdder.H"
|
#include "fvMeshAdder.H"
|
||||||
#include "polyTopoChange.H"
|
#include "polyTopoChange.H"
|
||||||
#include "extrapolatedCalculatedFvPatchFields.H"
|
#include "extrapolatedCalculatedFvPatchFields.H"
|
||||||
|
#include "regionProperties.H"
|
||||||
|
|
||||||
using namespace Foam;
|
using namespace Foam;
|
||||||
|
|
||||||
@ -479,6 +480,7 @@ int main(int argc, char *argv[])
|
|||||||
);
|
);
|
||||||
|
|
||||||
#include "addRegionOption.H"
|
#include "addRegionOption.H"
|
||||||
|
#include "addAllRegionsOption.H"
|
||||||
#include "setRootCase.H"
|
#include "setRootCase.H"
|
||||||
#include "createTime.H"
|
#include "createTime.H"
|
||||||
|
|
||||||
@ -498,18 +500,21 @@ int main(int argc, char *argv[])
|
|||||||
<< endl;
|
<< endl;
|
||||||
|
|
||||||
|
|
||||||
word regionName = polyMesh::defaultRegion;
|
const wordList regionNames(selectRegionNames(args, runTime));
|
||||||
word regionDir = word::null;
|
if (regionNames.size() > 1)
|
||||||
|
|
||||||
if
|
|
||||||
(
|
|
||||||
args.optionReadIfPresent("region", regionName)
|
|
||||||
&& regionName != polyMesh::defaultRegion
|
|
||||||
)
|
|
||||||
{
|
{
|
||||||
regionDir = regionName;
|
Info<< "Operating on regions " << regionNames[0];
|
||||||
Info<< "Operating on region " << regionName << nl << endl;
|
for (label regioni = 1; regioni < regionNames.size() - 1; ++ regioni)
|
||||||
|
{
|
||||||
|
Info<< ", " << regionNames[regioni];
|
||||||
|
}
|
||||||
|
Info<< " and " << regionNames.last() << nl << endl;
|
||||||
}
|
}
|
||||||
|
else if (regionNames[0] != polyMesh::defaultRegion)
|
||||||
|
{
|
||||||
|
Info<< "Operating on region " << regionNames[0] << nl << endl;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
scalar mergeTol = defaultMergeTol;
|
scalar mergeTol = defaultMergeTol;
|
||||||
args.optionReadIfPresent("mergeTol", mergeTol);
|
args.optionReadIfPresent("mergeTol", mergeTol);
|
||||||
@ -596,77 +601,288 @@ int main(int argc, char *argv[])
|
|||||||
databases[proci].setTime(timeDirs[timeI], timeI);
|
databases[proci].setTime(timeDirs[timeI], timeI);
|
||||||
}
|
}
|
||||||
|
|
||||||
IOobject facesIO
|
forAll(regionNames, regioni)
|
||||||
(
|
|
||||||
"faces",
|
|
||||||
databases[0].timeName(),
|
|
||||||
regionDir/polyMesh::meshSubDir,
|
|
||||||
databases[0],
|
|
||||||
IOobject::NO_READ,
|
|
||||||
IOobject::NO_WRITE
|
|
||||||
);
|
|
||||||
|
|
||||||
|
|
||||||
// Problem: faceCompactIOList recognises both 'faceList' and
|
|
||||||
// 'faceCompactList' so we should be lenient when doing
|
|
||||||
// typeHeaderOk
|
|
||||||
if (!facesIO.typeHeaderOk<faceCompactIOList>(false))
|
|
||||||
{
|
{
|
||||||
Info<< "No mesh." << nl << endl;
|
const word& regionName = regionNames[regioni];
|
||||||
continue;
|
const word regionDir =
|
||||||
}
|
regionName == polyMesh::defaultRegion
|
||||||
|
? word::null
|
||||||
|
: regionName;
|
||||||
|
|
||||||
|
IOobject facesIO
|
||||||
|
(
|
||||||
|
"faces",
|
||||||
|
databases[0].timeName(),
|
||||||
|
regionDir/polyMesh::meshSubDir,
|
||||||
|
databases[0],
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::NO_WRITE
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
// Read point on individual processors to determine merge tolerance
|
// Problem: faceCompactIOList recognises both 'faceList' and
|
||||||
// (otherwise single cell domains might give problems)
|
// 'faceCompactList' so we should be lenient when doing
|
||||||
|
// typeHeaderOk
|
||||||
const boundBox bb = procBounds(args, databases, regionDir);
|
if (!facesIO.typeHeaderOk<faceCompactIOList>(false))
|
||||||
const scalar mergeDist = mergeTol*bb.mag();
|
|
||||||
|
|
||||||
Info<< "Overall mesh bounding box : " << bb << nl
|
|
||||||
<< "Relative tolerance : " << mergeTol << nl
|
|
||||||
<< "Absolute matching distance : " << mergeDist << nl
|
|
||||||
<< endl;
|
|
||||||
|
|
||||||
|
|
||||||
// Addressing from processor to reconstructed case
|
|
||||||
labelListList cellProcAddressing(nProcs);
|
|
||||||
labelListList faceProcAddressing(nProcs);
|
|
||||||
labelListList pointProcAddressing(nProcs);
|
|
||||||
labelListList boundaryProcAddressing(nProcs);
|
|
||||||
|
|
||||||
// Internal faces on the final reconstructed mesh
|
|
||||||
label masterInternalFaces;
|
|
||||||
|
|
||||||
// Owner addressing on the final reconstructed mesh
|
|
||||||
labelList masterOwner;
|
|
||||||
|
|
||||||
{
|
|
||||||
// Construct empty mesh.
|
|
||||||
// fvMesh** masterMesh = new fvMesh*[nProcs];
|
|
||||||
PtrList<fvMesh> masterMesh(nProcs);
|
|
||||||
|
|
||||||
for (label proci=0; proci<nProcs; proci++)
|
|
||||||
{
|
{
|
||||||
masterMesh.set
|
Info<< "No mesh." << nl << endl;
|
||||||
(
|
continue;
|
||||||
proci,
|
}
|
||||||
new fvMesh
|
|
||||||
|
|
||||||
|
// Read point on individual processors to determine merge tolerance
|
||||||
|
// (otherwise single cell domains might give problems)
|
||||||
|
|
||||||
|
const boundBox bb = procBounds(args, databases, regionDir);
|
||||||
|
const scalar mergeDist = mergeTol*bb.mag();
|
||||||
|
|
||||||
|
Info<< "Overall mesh bounding box : " << bb << nl
|
||||||
|
<< "Relative tolerance : " << mergeTol << nl
|
||||||
|
<< "Absolute matching distance : " << mergeDist << nl
|
||||||
|
<< endl;
|
||||||
|
|
||||||
|
|
||||||
|
// Addressing from processor to reconstructed case
|
||||||
|
labelListList cellProcAddressing(nProcs);
|
||||||
|
labelListList faceProcAddressing(nProcs);
|
||||||
|
labelListList pointProcAddressing(nProcs);
|
||||||
|
labelListList boundaryProcAddressing(nProcs);
|
||||||
|
|
||||||
|
// Internal faces on the final reconstructed mesh
|
||||||
|
label masterInternalFaces;
|
||||||
|
|
||||||
|
// Owner addressing on the final reconstructed mesh
|
||||||
|
labelList masterOwner;
|
||||||
|
|
||||||
|
{
|
||||||
|
// Construct empty mesh.
|
||||||
|
// fvMesh** masterMesh = new fvMesh*[nProcs];
|
||||||
|
PtrList<fvMesh> masterMesh(nProcs);
|
||||||
|
|
||||||
|
for (label proci=0; proci<nProcs; proci++)
|
||||||
|
{
|
||||||
|
masterMesh.set
|
||||||
|
(
|
||||||
|
proci,
|
||||||
|
new fvMesh
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
regionName,
|
||||||
|
runTime.timeName(),
|
||||||
|
runTime,
|
||||||
|
IOobject::NO_READ
|
||||||
|
),
|
||||||
|
xferCopy(pointField()),
|
||||||
|
xferCopy(faceList()),
|
||||||
|
xferCopy(cellList())
|
||||||
|
)
|
||||||
|
);
|
||||||
|
|
||||||
|
fvMesh meshToAdd
|
||||||
(
|
(
|
||||||
IOobject
|
IOobject
|
||||||
(
|
(
|
||||||
regionName,
|
regionName,
|
||||||
runTime.timeName(),
|
databases[proci].timeName(),
|
||||||
runTime,
|
databases[proci]
|
||||||
IOobject::NO_READ
|
)
|
||||||
),
|
);
|
||||||
xferCopy(pointField()),
|
|
||||||
xferCopy(faceList()),
|
// Initialize its addressing
|
||||||
xferCopy(cellList())
|
cellProcAddressing[proci] = identity(meshToAdd.nCells());
|
||||||
)
|
faceProcAddressing[proci] = identity(meshToAdd.nFaces());
|
||||||
|
pointProcAddressing[proci] = identity(meshToAdd.nPoints());
|
||||||
|
boundaryProcAddressing[proci] =
|
||||||
|
identity(meshToAdd.boundaryMesh().size());
|
||||||
|
|
||||||
|
// Find geometrically shared points/faces.
|
||||||
|
autoPtr<faceCoupleInfo> couples = determineCoupledFaces
|
||||||
|
(
|
||||||
|
fullMatch,
|
||||||
|
proci,
|
||||||
|
proci,
|
||||||
|
masterMesh[proci],
|
||||||
|
proci,
|
||||||
|
proci,
|
||||||
|
meshToAdd,
|
||||||
|
mergeDist
|
||||||
|
);
|
||||||
|
|
||||||
|
// Add elements to mesh
|
||||||
|
autoPtr<mapAddedPolyMesh> map = fvMeshAdder::add
|
||||||
|
(
|
||||||
|
masterMesh[proci],
|
||||||
|
meshToAdd,
|
||||||
|
couples
|
||||||
|
);
|
||||||
|
|
||||||
|
// Added processor
|
||||||
|
renumber(map().addedCellMap(), cellProcAddressing[proci]);
|
||||||
|
renumber(map().addedFaceMap(), faceProcAddressing[proci]);
|
||||||
|
renumber(map().addedPointMap(), pointProcAddressing[proci]);
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().addedPatchMap(),
|
||||||
|
boundaryProcAddressing[proci]
|
||||||
|
);
|
||||||
|
}
|
||||||
|
for (label step=2; step<nProcs*2; step*=2)
|
||||||
|
{
|
||||||
|
for (label proci=0; proci<nProcs; proci+=step)
|
||||||
|
{
|
||||||
|
label next = proci + step/2;
|
||||||
|
if(next >= nProcs)
|
||||||
|
{
|
||||||
|
continue;
|
||||||
|
}
|
||||||
|
|
||||||
|
Info<< "Merging mesh " << proci << " with " << next
|
||||||
|
<< endl;
|
||||||
|
|
||||||
|
// Find geometrically shared points/faces.
|
||||||
|
autoPtr<faceCoupleInfo> couples = determineCoupledFaces
|
||||||
|
(
|
||||||
|
fullMatch,
|
||||||
|
proci,
|
||||||
|
next,
|
||||||
|
masterMesh[proci],
|
||||||
|
next,
|
||||||
|
proci+step,
|
||||||
|
masterMesh[next],
|
||||||
|
mergeDist
|
||||||
|
);
|
||||||
|
|
||||||
|
// Add elements to mesh
|
||||||
|
autoPtr<mapAddedPolyMesh> map = fvMeshAdder::add
|
||||||
|
(
|
||||||
|
masterMesh[proci],
|
||||||
|
masterMesh[next],
|
||||||
|
couples
|
||||||
|
);
|
||||||
|
|
||||||
|
// Processors that were already in masterMesh
|
||||||
|
for (label mergedI=proci; mergedI<next; mergedI++)
|
||||||
|
{
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().oldCellMap(),
|
||||||
|
cellProcAddressing[mergedI]
|
||||||
|
);
|
||||||
|
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().oldFaceMap(),
|
||||||
|
faceProcAddressing[mergedI]
|
||||||
|
);
|
||||||
|
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().oldPointMap(),
|
||||||
|
pointProcAddressing[mergedI]
|
||||||
|
);
|
||||||
|
|
||||||
|
// Note: boundary is special since can contain -1.
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().oldPatchMap(),
|
||||||
|
boundaryProcAddressing[mergedI]
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
// Added processor
|
||||||
|
for
|
||||||
|
(
|
||||||
|
label addedI=next;
|
||||||
|
addedI<min(proci+step, nProcs);
|
||||||
|
addedI++
|
||||||
|
)
|
||||||
|
{
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().addedCellMap(),
|
||||||
|
cellProcAddressing[addedI]
|
||||||
|
);
|
||||||
|
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().addedFaceMap(),
|
||||||
|
faceProcAddressing[addedI]
|
||||||
|
);
|
||||||
|
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().addedPointMap(),
|
||||||
|
pointProcAddressing[addedI]
|
||||||
|
);
|
||||||
|
|
||||||
|
renumber
|
||||||
|
(
|
||||||
|
map().addedPatchMap(),
|
||||||
|
boundaryProcAddressing[addedI]
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
masterMesh.set(next, nullptr);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
for (label proci=0; proci<nProcs; proci++)
|
||||||
|
{
|
||||||
|
Info<< "Reading mesh to add from "
|
||||||
|
<< databases[proci].caseName()
|
||||||
|
<< " for time = " << databases[proci].timeName()
|
||||||
|
<< nl << nl << endl;
|
||||||
|
}
|
||||||
|
|
||||||
|
// See if any points on the mastermesh have become connected
|
||||||
|
// because of connections through processor meshes.
|
||||||
|
mergeSharedPoints
|
||||||
|
(
|
||||||
|
mergeDist,
|
||||||
|
masterMesh[0],
|
||||||
|
pointProcAddressing
|
||||||
);
|
);
|
||||||
|
|
||||||
fvMesh meshToAdd
|
// Save some properties on the reconstructed mesh
|
||||||
|
masterInternalFaces = masterMesh[0].nInternalFaces();
|
||||||
|
masterOwner = masterMesh[0].faceOwner();
|
||||||
|
|
||||||
|
|
||||||
|
Info<< "\nWriting merged mesh to "
|
||||||
|
<< runTime.path()/runTime.timeName()
|
||||||
|
<< nl << endl;
|
||||||
|
|
||||||
|
if (!masterMesh[0].write())
|
||||||
|
{
|
||||||
|
FatalErrorInFunction
|
||||||
|
<< "Failed writing polyMesh."
|
||||||
|
<< exit(FatalError);
|
||||||
|
}
|
||||||
|
|
||||||
|
if (writeCellDist)
|
||||||
|
{
|
||||||
|
writeCellDistance
|
||||||
|
(
|
||||||
|
runTime,
|
||||||
|
masterMesh[0],
|
||||||
|
cellProcAddressing
|
||||||
|
);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// Write the addressing
|
||||||
|
|
||||||
|
Info<< "Reconstructing the addressing from the processor meshes"
|
||||||
|
<< " to the newly reconstructed mesh" << nl << endl;
|
||||||
|
|
||||||
|
forAll(databases, proci)
|
||||||
|
{
|
||||||
|
Info<< "Reading processor " << proci << " mesh from "
|
||||||
|
<< databases[proci].caseName() << endl;
|
||||||
|
|
||||||
|
polyMesh procMesh
|
||||||
(
|
(
|
||||||
IOobject
|
IOobject
|
||||||
(
|
(
|
||||||
@ -676,336 +892,144 @@ int main(int argc, char *argv[])
|
|||||||
)
|
)
|
||||||
);
|
);
|
||||||
|
|
||||||
// Initialize its addressing
|
|
||||||
cellProcAddressing[proci] = identity(meshToAdd.nCells());
|
|
||||||
faceProcAddressing[proci] = identity(meshToAdd.nFaces());
|
|
||||||
pointProcAddressing[proci] = identity(meshToAdd.nPoints());
|
|
||||||
boundaryProcAddressing[proci] =
|
|
||||||
identity(meshToAdd.boundaryMesh().size());
|
|
||||||
|
|
||||||
// Find geometrically shared points/faces.
|
// From processor point to reconstructed mesh point
|
||||||
autoPtr<faceCoupleInfo> couples = determineCoupledFaces
|
|
||||||
|
Info<< "Writing pointProcAddressing to "
|
||||||
|
<< databases[proci].caseName()
|
||||||
|
/procMesh.facesInstance()
|
||||||
|
/polyMesh::meshSubDir
|
||||||
|
<< endl;
|
||||||
|
|
||||||
|
labelIOList
|
||||||
(
|
(
|
||||||
fullMatch,
|
IOobject
|
||||||
proci,
|
(
|
||||||
proci,
|
"pointProcAddressing",
|
||||||
masterMesh[proci],
|
procMesh.facesInstance(),
|
||||||
proci,
|
polyMesh::meshSubDir,
|
||||||
proci,
|
procMesh,
|
||||||
meshToAdd,
|
IOobject::NO_READ,
|
||||||
mergeDist
|
IOobject::NO_WRITE,
|
||||||
|
false // Do not register
|
||||||
|
),
|
||||||
|
pointProcAddressing[proci]
|
||||||
|
).write();
|
||||||
|
|
||||||
|
|
||||||
|
// From processor face to reconstructed mesh face
|
||||||
|
|
||||||
|
Info<< "Writing faceProcAddressing to "
|
||||||
|
<< databases[proci].caseName()
|
||||||
|
/procMesh.facesInstance()
|
||||||
|
/polyMesh::meshSubDir
|
||||||
|
<< endl;
|
||||||
|
|
||||||
|
labelIOList faceProcAddr
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
"faceProcAddressing",
|
||||||
|
procMesh.facesInstance(),
|
||||||
|
polyMesh::meshSubDir,
|
||||||
|
procMesh,
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::NO_WRITE,
|
||||||
|
false // Do not register
|
||||||
|
),
|
||||||
|
faceProcAddressing[proci]
|
||||||
);
|
);
|
||||||
|
|
||||||
// Add elements to mesh
|
// Now add turning index to faceProcAddressing.
|
||||||
autoPtr<mapAddedPolyMesh> map = fvMeshAdder::add
|
// See reconstructPar for meaning of turning index.
|
||||||
(
|
forAll(faceProcAddr, procFacei)
|
||||||
masterMesh[proci],
|
|
||||||
meshToAdd,
|
|
||||||
couples
|
|
||||||
);
|
|
||||||
|
|
||||||
// Added processor
|
|
||||||
renumber(map().addedCellMap(), cellProcAddressing[proci]);
|
|
||||||
renumber(map().addedFaceMap(), faceProcAddressing[proci]);
|
|
||||||
renumber(map().addedPointMap(), pointProcAddressing[proci]);
|
|
||||||
renumber(map().addedPatchMap(), boundaryProcAddressing[proci]);
|
|
||||||
}
|
|
||||||
for (label step=2; step<nProcs*2; step*=2)
|
|
||||||
{
|
|
||||||
for (label proci=0; proci<nProcs; proci+=step)
|
|
||||||
{
|
{
|
||||||
label next = proci + step/2;
|
const label masterFacei = faceProcAddr[procFacei];
|
||||||
if(next >= nProcs)
|
|
||||||
{
|
|
||||||
continue;
|
|
||||||
}
|
|
||||||
|
|
||||||
Info<< "Merging mesh " << proci << " with " << next << endl;
|
if
|
||||||
|
|
||||||
// Find geometrically shared points/faces.
|
|
||||||
autoPtr<faceCoupleInfo> couples = determineCoupledFaces
|
|
||||||
(
|
(
|
||||||
fullMatch,
|
!procMesh.isInternalFace(procFacei)
|
||||||
proci,
|
&& masterFacei < masterInternalFaces
|
||||||
next,
|
|
||||||
masterMesh[proci],
|
|
||||||
next,
|
|
||||||
proci+step,
|
|
||||||
masterMesh[next],
|
|
||||||
mergeDist
|
|
||||||
);
|
|
||||||
|
|
||||||
// Add elements to mesh
|
|
||||||
autoPtr<mapAddedPolyMesh> map = fvMeshAdder::add
|
|
||||||
(
|
|
||||||
masterMesh[proci],
|
|
||||||
masterMesh[next],
|
|
||||||
couples
|
|
||||||
);
|
|
||||||
|
|
||||||
// Processors that were already in masterMesh
|
|
||||||
for (label mergedI=proci; mergedI<next; mergedI++)
|
|
||||||
{
|
|
||||||
renumber
|
|
||||||
(
|
|
||||||
map().oldCellMap(),
|
|
||||||
cellProcAddressing[mergedI]
|
|
||||||
);
|
|
||||||
|
|
||||||
renumber
|
|
||||||
(
|
|
||||||
map().oldFaceMap(),
|
|
||||||
faceProcAddressing[mergedI]
|
|
||||||
);
|
|
||||||
|
|
||||||
renumber
|
|
||||||
(
|
|
||||||
map().oldPointMap(),
|
|
||||||
pointProcAddressing[mergedI]
|
|
||||||
);
|
|
||||||
|
|
||||||
// Note: boundary is special since can contain -1.
|
|
||||||
renumber
|
|
||||||
(
|
|
||||||
map().oldPatchMap(),
|
|
||||||
boundaryProcAddressing[mergedI]
|
|
||||||
);
|
|
||||||
}
|
|
||||||
|
|
||||||
// Added processor
|
|
||||||
for
|
|
||||||
(
|
|
||||||
label addedI=next;
|
|
||||||
addedI<min(proci+step, nProcs);
|
|
||||||
addedI++
|
|
||||||
)
|
)
|
||||||
{
|
{
|
||||||
renumber
|
// proc face is now external but used to be internal
|
||||||
(
|
// face. Check if we have owner or neighbour.
|
||||||
map().addedCellMap(),
|
|
||||||
cellProcAddressing[addedI]
|
|
||||||
);
|
|
||||||
|
|
||||||
renumber
|
label procOwn = procMesh.faceOwner()[procFacei];
|
||||||
(
|
label masterOwn = masterOwner[masterFacei];
|
||||||
map().addedFaceMap(),
|
|
||||||
faceProcAddressing[addedI]
|
|
||||||
);
|
|
||||||
|
|
||||||
renumber
|
if (cellProcAddressing[proci][procOwn] == masterOwn)
|
||||||
(
|
{
|
||||||
map().addedPointMap(),
|
// No turning. Offset by 1.
|
||||||
pointProcAddressing[addedI]
|
faceProcAddr[procFacei]++;
|
||||||
);
|
}
|
||||||
|
else
|
||||||
renumber
|
{
|
||||||
(
|
// Turned face.
|
||||||
map().addedPatchMap(),
|
faceProcAddr[procFacei] =
|
||||||
boundaryProcAddressing[addedI]
|
-1 - faceProcAddr[procFacei];
|
||||||
);
|
}
|
||||||
}
|
}
|
||||||
|
else
|
||||||
masterMesh.set(next, nullptr);
|
|
||||||
}
|
|
||||||
}
|
|
||||||
|
|
||||||
for (label proci=0; proci<nProcs; proci++)
|
|
||||||
{
|
|
||||||
Info<< "Reading mesh to add from "
|
|
||||||
<< databases[proci].caseName()
|
|
||||||
<< " for time = " << databases[proci].timeName()
|
|
||||||
<< nl << nl << endl;
|
|
||||||
}
|
|
||||||
|
|
||||||
// See if any points on the mastermesh have become connected
|
|
||||||
// because of connections through processor meshes.
|
|
||||||
mergeSharedPoints(mergeDist, masterMesh[0], pointProcAddressing);
|
|
||||||
|
|
||||||
// Save some properties on the reconstructed mesh
|
|
||||||
masterInternalFaces = masterMesh[0].nInternalFaces();
|
|
||||||
masterOwner = masterMesh[0].faceOwner();
|
|
||||||
|
|
||||||
|
|
||||||
Info<< "\nWriting merged mesh to "
|
|
||||||
<< runTime.path()/runTime.timeName()
|
|
||||||
<< nl << endl;
|
|
||||||
|
|
||||||
if (!masterMesh[0].write())
|
|
||||||
{
|
|
||||||
FatalErrorInFunction
|
|
||||||
<< "Failed writing polyMesh."
|
|
||||||
<< exit(FatalError);
|
|
||||||
}
|
|
||||||
|
|
||||||
if (writeCellDist)
|
|
||||||
{
|
|
||||||
writeCellDistance
|
|
||||||
(
|
|
||||||
runTime,
|
|
||||||
masterMesh[0],
|
|
||||||
cellProcAddressing
|
|
||||||
);
|
|
||||||
}
|
|
||||||
}
|
|
||||||
|
|
||||||
|
|
||||||
// Write the addressing
|
|
||||||
|
|
||||||
Info<< "Reconstructing the addressing from the processor meshes"
|
|
||||||
<< " to the newly reconstructed mesh" << nl << endl;
|
|
||||||
|
|
||||||
forAll(databases, proci)
|
|
||||||
{
|
|
||||||
Info<< "Reading processor " << proci << " mesh from "
|
|
||||||
<< databases[proci].caseName() << endl;
|
|
||||||
|
|
||||||
polyMesh procMesh
|
|
||||||
(
|
|
||||||
IOobject
|
|
||||||
(
|
|
||||||
regionName,
|
|
||||||
databases[proci].timeName(),
|
|
||||||
databases[proci]
|
|
||||||
)
|
|
||||||
);
|
|
||||||
|
|
||||||
|
|
||||||
// From processor point to reconstructed mesh point
|
|
||||||
|
|
||||||
Info<< "Writing pointProcAddressing to "
|
|
||||||
<< databases[proci].caseName()
|
|
||||||
/procMesh.facesInstance()
|
|
||||||
/polyMesh::meshSubDir
|
|
||||||
<< endl;
|
|
||||||
|
|
||||||
labelIOList
|
|
||||||
(
|
|
||||||
IOobject
|
|
||||||
(
|
|
||||||
"pointProcAddressing",
|
|
||||||
procMesh.facesInstance(),
|
|
||||||
polyMesh::meshSubDir,
|
|
||||||
procMesh,
|
|
||||||
IOobject::NO_READ,
|
|
||||||
IOobject::NO_WRITE,
|
|
||||||
false // Do not register
|
|
||||||
),
|
|
||||||
pointProcAddressing[proci]
|
|
||||||
).write();
|
|
||||||
|
|
||||||
|
|
||||||
// From processor face to reconstructed mesh face
|
|
||||||
|
|
||||||
Info<< "Writing faceProcAddressing to "
|
|
||||||
<< databases[proci].caseName()
|
|
||||||
/procMesh.facesInstance()
|
|
||||||
/polyMesh::meshSubDir
|
|
||||||
<< endl;
|
|
||||||
|
|
||||||
labelIOList faceProcAddr
|
|
||||||
(
|
|
||||||
IOobject
|
|
||||||
(
|
|
||||||
"faceProcAddressing",
|
|
||||||
procMesh.facesInstance(),
|
|
||||||
polyMesh::meshSubDir,
|
|
||||||
procMesh,
|
|
||||||
IOobject::NO_READ,
|
|
||||||
IOobject::NO_WRITE,
|
|
||||||
false // Do not register
|
|
||||||
),
|
|
||||||
faceProcAddressing[proci]
|
|
||||||
);
|
|
||||||
|
|
||||||
// Now add turning index to faceProcAddressing.
|
|
||||||
// See reconstructPar for meaning of turning index.
|
|
||||||
forAll(faceProcAddr, procFacei)
|
|
||||||
{
|
|
||||||
const label masterFacei = faceProcAddr[procFacei];
|
|
||||||
|
|
||||||
if
|
|
||||||
(
|
|
||||||
!procMesh.isInternalFace(procFacei)
|
|
||||||
&& masterFacei < masterInternalFaces
|
|
||||||
)
|
|
||||||
{
|
|
||||||
// proc face is now external but used to be internal face.
|
|
||||||
// Check if we have owner or neighbour.
|
|
||||||
|
|
||||||
label procOwn = procMesh.faceOwner()[procFacei];
|
|
||||||
label masterOwn = masterOwner[masterFacei];
|
|
||||||
|
|
||||||
if (cellProcAddressing[proci][procOwn] == masterOwn)
|
|
||||||
{
|
{
|
||||||
// No turning. Offset by 1.
|
// No turning. Offset by 1.
|
||||||
faceProcAddr[procFacei]++;
|
faceProcAddr[procFacei]++;
|
||||||
}
|
}
|
||||||
else
|
|
||||||
{
|
|
||||||
// Turned face.
|
|
||||||
faceProcAddr[procFacei] =
|
|
||||||
-1 - faceProcAddr[procFacei];
|
|
||||||
}
|
|
||||||
}
|
|
||||||
else
|
|
||||||
{
|
|
||||||
// No turning. Offset by 1.
|
|
||||||
faceProcAddr[procFacei]++;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
faceProcAddr.write();
|
||||||
|
|
||||||
|
|
||||||
|
// From processor cell to reconstructed mesh cell
|
||||||
|
|
||||||
|
Info<< "Writing cellProcAddressing to "
|
||||||
|
<< databases[proci].caseName()
|
||||||
|
/procMesh.facesInstance()
|
||||||
|
/polyMesh::meshSubDir
|
||||||
|
<< endl;
|
||||||
|
|
||||||
|
labelIOList
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
"cellProcAddressing",
|
||||||
|
procMesh.facesInstance(),
|
||||||
|
polyMesh::meshSubDir,
|
||||||
|
procMesh,
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::NO_WRITE,
|
||||||
|
false // Do not register
|
||||||
|
),
|
||||||
|
cellProcAddressing[proci]
|
||||||
|
).write();
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
// From processor patch to reconstructed mesh patch
|
||||||
|
|
||||||
|
Info<< "Writing boundaryProcAddressing to "
|
||||||
|
<< databases[proci].caseName()
|
||||||
|
/procMesh.facesInstance()
|
||||||
|
/polyMesh::meshSubDir
|
||||||
|
<< endl;
|
||||||
|
|
||||||
|
labelIOList
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
"boundaryProcAddressing",
|
||||||
|
procMesh.facesInstance(),
|
||||||
|
polyMesh::meshSubDir,
|
||||||
|
procMesh,
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::NO_WRITE,
|
||||||
|
false // Do not register
|
||||||
|
),
|
||||||
|
boundaryProcAddressing[proci]
|
||||||
|
).write();
|
||||||
|
|
||||||
|
Info<< endl;
|
||||||
}
|
}
|
||||||
|
|
||||||
faceProcAddr.write();
|
|
||||||
|
|
||||||
|
|
||||||
// From processor cell to reconstructed mesh cell
|
|
||||||
|
|
||||||
Info<< "Writing cellProcAddressing to "
|
|
||||||
<< databases[proci].caseName()
|
|
||||||
/procMesh.facesInstance()
|
|
||||||
/polyMesh::meshSubDir
|
|
||||||
<< endl;
|
|
||||||
|
|
||||||
labelIOList
|
|
||||||
(
|
|
||||||
IOobject
|
|
||||||
(
|
|
||||||
"cellProcAddressing",
|
|
||||||
procMesh.facesInstance(),
|
|
||||||
polyMesh::meshSubDir,
|
|
||||||
procMesh,
|
|
||||||
IOobject::NO_READ,
|
|
||||||
IOobject::NO_WRITE,
|
|
||||||
false // Do not register
|
|
||||||
),
|
|
||||||
cellProcAddressing[proci]
|
|
||||||
).write();
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
// From processor patch to reconstructed mesh patch
|
|
||||||
|
|
||||||
Info<< "Writing boundaryProcAddressing to "
|
|
||||||
<< databases[proci].caseName()
|
|
||||||
/procMesh.facesInstance()
|
|
||||||
/polyMesh::meshSubDir
|
|
||||||
<< endl;
|
|
||||||
|
|
||||||
labelIOList
|
|
||||||
(
|
|
||||||
IOobject
|
|
||||||
(
|
|
||||||
"boundaryProcAddressing",
|
|
||||||
procMesh.facesInstance(),
|
|
||||||
polyMesh::meshSubDir,
|
|
||||||
procMesh,
|
|
||||||
IOobject::NO_READ,
|
|
||||||
IOobject::NO_WRITE,
|
|
||||||
false // Do not register
|
|
||||||
),
|
|
||||||
boundaryProcAddressing[proci]
|
|
||||||
).write();
|
|
||||||
|
|
||||||
Info<< endl;
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|||||||
@ -6,6 +6,8 @@ cd ${0%/*} || exit 1 # Run from this directory
|
|||||||
|
|
||||||
cleanCase
|
cleanCase
|
||||||
rm -rf \
|
rm -rf \
|
||||||
|
0/cellToRegion \
|
||||||
|
0/*/cellToRegion \
|
||||||
constant/*/polyMesh \
|
constant/*/polyMesh \
|
||||||
constant/extendedFeatureEdgeMesh \
|
constant/extendedFeatureEdgeMesh \
|
||||||
constant/triSurface/*.eMesh
|
constant/triSurface/*.eMesh
|
||||||
|
|||||||
@ -7,17 +7,13 @@ cd ${0%/*} || exit 1 # Run from this directory
|
|||||||
# Meshing
|
# Meshing
|
||||||
runApplication blockMesh
|
runApplication blockMesh
|
||||||
runApplication topoSet
|
runApplication topoSet
|
||||||
runApplication decomposePar -copyZero
|
runApplication splitMeshRegions -cellZones -overwrite
|
||||||
runParallel splitMeshRegions -cellZones -overwrite
|
runApplication decomposePar -copyZero -allRegions
|
||||||
|
|
||||||
# Simulation
|
# Simulation
|
||||||
runParallel $(getApplication)
|
runParallel $(getApplication)
|
||||||
|
|
||||||
# Reconstruct
|
# Reconstruct
|
||||||
for region in bottomWater topAir heater leftSolid rightSolid
|
|
||||||
do
|
|
||||||
runApplication -s $region reconstructParMesh -constant -region $region
|
|
||||||
done
|
|
||||||
runApplication reconstructPar -allRegions
|
runApplication reconstructPar -allRegions
|
||||||
|
|
||||||
# Post-process
|
# Post-process
|
||||||
|
|||||||
@ -15,10 +15,7 @@ runParallel splitMeshRegions -cellZones -overwrite
|
|||||||
runParallel $(getApplication)
|
runParallel $(getApplication)
|
||||||
|
|
||||||
# Reconstruct
|
# Reconstruct
|
||||||
for region in bottomAir topAir heater leftSolid rightSolid
|
runApplication reconstructParMesh -constant -allRegions
|
||||||
do
|
|
||||||
runApplication -s $region reconstructParMesh -constant -region $region
|
|
||||||
done
|
|
||||||
runApplication reconstructPar -allRegions
|
runApplication reconstructPar -allRegions
|
||||||
|
|
||||||
# Post-process
|
# Post-process
|
||||||
|
|||||||
Reference in New Issue
Block a user