Resolve cell and point level load balancing for cases without adaptive refinement
This commit is contained in:
parent
d23e612f80
commit
8a3e798d5a
1 changed files with 134 additions and 84 deletions
|
@ -164,13 +164,13 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
|
||||
// Distribute cell and point level for AMR + DLB runs. VV, 18/May/2018
|
||||
|
||||
// Get (non-const) reference to cellLevel
|
||||
labelIOField& cellLevel = const_cast<labelIOField&>
|
||||
(this->lookupObject<labelIOField>("cellLevel"));
|
||||
|
||||
// Get (non-const) reference to pointLevel
|
||||
labelIOField& pointLevel = const_cast<labelIOField&>
|
||||
(this->lookupObject<labelIOField>("pointLevel"));
|
||||
// Check cellLevel and pointLevel
|
||||
// Note: possible sync problem if cellLevel or pointLevel is not present
|
||||
// on all processors in comms. However, parallel reduce check
|
||||
// is not allowed
|
||||
// HJ, 22/Oct/2018
|
||||
bool cellLevelFound = this->foundObject<labelIOField>("cellLevel");
|
||||
bool pointLevelFound = this->foundObject<labelIOField>("pointLevel");
|
||||
|
||||
//HJ, HERE: remove the fields that should not be load balanced
|
||||
|
||||
|
@ -245,7 +245,7 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
for (label procI = 0; procI < meshDecomp.nProcs(); procI++)
|
||||
{
|
||||
// Check if there is a mesh to send
|
||||
if (curMigratedCells[procI] > 0)
|
||||
if (migratedCells[Pstream::myProcNo()][procI] > 0)
|
||||
{
|
||||
// Send mesh and field to processor procI
|
||||
|
||||
|
@ -273,6 +273,11 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
{
|
||||
// Pout<< "Send mesh and fields to processor " << procI << endl;
|
||||
|
||||
// OFstream toProc
|
||||
// (
|
||||
// "from" + Foam::name(Pstream::myProcNo())
|
||||
// + "To" + Foam::name(procI)
|
||||
// );
|
||||
OPstream toProc
|
||||
(
|
||||
Pstream::blocking,
|
||||
|
@ -289,6 +294,11 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
sendFields(surfaceVectorFields, fieldDecomposer, toProc);
|
||||
|
||||
// Send cell level with procCellAddressing
|
||||
if (cellLevelFound)
|
||||
{
|
||||
const labelIOField& cellLevel =
|
||||
this->lookupObject<labelIOField>("cellLevel");
|
||||
|
||||
toProc <<
|
||||
labelList
|
||||
(
|
||||
|
@ -296,6 +306,12 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
meshDecomp.procCellAddressing()[procI]
|
||||
)
|
||||
<< nl;
|
||||
}
|
||||
|
||||
if (pointLevelFound)
|
||||
{
|
||||
const labelIOField& pointLevel =
|
||||
this->lookupObject<labelIOField>("pointLevel");
|
||||
|
||||
// Send point level with procPointAddressing
|
||||
toProc <<
|
||||
|
@ -305,6 +321,7 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
meshDecomp.procPointAddressing()[procI]
|
||||
)
|
||||
<< nl;
|
||||
}
|
||||
|
||||
// Send clouds
|
||||
forAll(cloudDistributes, cloudI)
|
||||
|
@ -355,6 +372,12 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
receivedSurfaceVectorFields
|
||||
);
|
||||
|
||||
// Insert cell level
|
||||
if (cellLevelFound)
|
||||
{
|
||||
const labelIOField& cellLevel =
|
||||
this->lookupObject<labelIOField>("cellLevel");
|
||||
|
||||
// Insert cell level
|
||||
receivedCellLevel.set
|
||||
(
|
||||
|
@ -366,6 +389,13 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
[Pstream::myProcNo()]
|
||||
)
|
||||
);
|
||||
}
|
||||
|
||||
// Insert point level
|
||||
if (pointLevelFound)
|
||||
{
|
||||
const labelIOField& pointLevel =
|
||||
this->lookupObject<labelIOField>("pointLevel");
|
||||
|
||||
// Insert point level
|
||||
receivedPointLevel.set
|
||||
|
@ -378,13 +408,13 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
[Pstream::myProcNo()]
|
||||
)
|
||||
);
|
||||
}
|
||||
|
||||
//HJ Insert clouds missing. HJ, 12/Oct/2018
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// Collect pieces of mesh and fields from other processors
|
||||
for (label procI = 0; procI < meshDecomp.nProcs(); procI++)
|
||||
{
|
||||
|
@ -393,10 +423,14 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
// Check if there is a mesh to send
|
||||
if (migratedCells[procI][Pstream::myProcNo()] > 0)
|
||||
{
|
||||
// Pout<< "Receive mesh and fields from " << procI
|
||||
// << endl;
|
||||
// Pout<< "Receive mesh and fields from " << procI << endl;
|
||||
|
||||
// Note: communication can be optimised. HJ, 27/Feb/2018
|
||||
// IFstream fromProc
|
||||
// (
|
||||
// "from" + Foam::name(procI)
|
||||
// + "To" + Foam::name(Pstream::myProcNo())
|
||||
// );
|
||||
IPstream fromProc
|
||||
(
|
||||
Pstream::blocking,
|
||||
|
@ -457,18 +491,24 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
);
|
||||
|
||||
// Receive cell level
|
||||
if (cellLevelFound)
|
||||
{
|
||||
receivedCellLevel.set
|
||||
(
|
||||
procI,
|
||||
new labelList(fromProc)
|
||||
);
|
||||
}
|
||||
|
||||
// Receive point level
|
||||
if (pointLevelFound)
|
||||
{
|
||||
receivedPointLevel.set
|
||||
(
|
||||
procI,
|
||||
new labelList(fromProc)
|
||||
);
|
||||
}
|
||||
|
||||
// Receive clouds
|
||||
forAll(cloudDistributes, cloudI)
|
||||
|
@ -710,6 +750,7 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
labelList patchStarts(patches.size());
|
||||
labelList oldPatchNMeshPoints(patches.size());
|
||||
labelListList patchPointMap(patches.size());
|
||||
|
||||
forAll(patches, patchI)
|
||||
{
|
||||
patchStarts[patchI] = patches[patchI].start();
|
||||
|
@ -757,7 +798,6 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
oldPatchNMeshPoints // oldPatchNMeshPoints
|
||||
);
|
||||
|
||||
|
||||
// Reset fvMesh and patches
|
||||
resetFvPrimitives
|
||||
(
|
||||
|
@ -816,10 +856,14 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
);
|
||||
|
||||
// Rebuild cell level field from components
|
||||
if (cellLevelFound)
|
||||
{
|
||||
// Get (non-const) reference to cellLevel
|
||||
labelIOField& cellLevel = const_cast<labelIOField&>
|
||||
(this->lookupObject<labelIOField>("cellLevel"));
|
||||
|
||||
if (cellLevel.size() != this->nCells())
|
||||
{
|
||||
Pout<< "Setting cell level size from: "
|
||||
<< cellLevel.size() << " to " << this->nCells() << endl;
|
||||
cellLevel.setSize(this->nCells());
|
||||
}
|
||||
|
||||
|
@ -834,12 +878,17 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
// Rebuild point level field from components
|
||||
if (pointLevelFound)
|
||||
{
|
||||
// Get (non-const) reference to pointLevel
|
||||
labelIOField& pointLevel = const_cast<labelIOField&>
|
||||
(this->lookupObject<labelIOField>("pointLevel"));
|
||||
|
||||
if (pointLevel.size() != this->nPoints())
|
||||
{
|
||||
Pout<< "Setting point level size from: "
|
||||
<< pointLevel.size() << " to " << this->nPoints() << endl;
|
||||
pointLevel.setSize(this->nPoints());
|
||||
}
|
||||
|
||||
|
@ -854,6 +903,7 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
// Rebuild clouds
|
||||
forAll(cloudDistributes, cloudI)
|
||||
|
@ -865,8 +915,8 @@ bool Foam::topoChangerFvMesh::loadBalance(const dictionary& decompDict)
|
|||
);
|
||||
}
|
||||
|
||||
// Debug
|
||||
checkMesh(true);
|
||||
// Debug: remove? HJ, 22/Oct/2018
|
||||
// checkMesh(true);
|
||||
|
||||
return true;
|
||||
}
|
||||
|
|
Reference in a new issue