Use labelField for refinement history to allow for parallel decomposition and reconstruction + load balancing
This commit is contained in:
parent
94af00fa04
commit
cffe702f75
3 changed files with 21 additions and 20 deletions
|
@ -1698,7 +1698,7 @@ Foam::hexRef8::hexRef8(const polyMesh& mesh)
|
|||
IOobject::READ_IF_PRESENT,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
labelList(mesh_.nCells(), 0)
|
||||
labelField(mesh_.nCells(), 0)
|
||||
),
|
||||
pointLevel_
|
||||
(
|
||||
|
@ -1711,7 +1711,7 @@ Foam::hexRef8::hexRef8(const polyMesh& mesh)
|
|||
IOobject::READ_IF_PRESENT,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
labelList(mesh_.nPoints(), 0)
|
||||
labelField(mesh_.nPoints(), 0)
|
||||
),
|
||||
level0Edge_(getLevel0EdgeLength()),
|
||||
history_
|
||||
|
@ -1779,8 +1779,8 @@ Foam::hexRef8::hexRef8(const polyMesh& mesh)
|
|||
Foam::hexRef8::hexRef8
|
||||
(
|
||||
const polyMesh& mesh,
|
||||
const labelList& cellLevel,
|
||||
const labelList& pointLevel,
|
||||
const labelField& cellLevel,
|
||||
const labelField& pointLevel,
|
||||
const refinementHistory& history
|
||||
)
|
||||
:
|
||||
|
@ -1833,8 +1833,8 @@ Foam::hexRef8::hexRef8
|
|||
{
|
||||
FatalErrorIn
|
||||
(
|
||||
"hexRef8::hexRef8(const polyMesh&, const labelList&"
|
||||
", const labelList&, const refinementHistory&)"
|
||||
"hexRef8::hexRef8(const polyMesh&, const labelField&"
|
||||
", const labelField&, const refinementHistory&)"
|
||||
) << "History enabled but number of visible cells in it "
|
||||
<< history_.visibleCells().size()
|
||||
<< " is not equal to the number of cells in the mesh "
|
||||
|
@ -1849,8 +1849,8 @@ Foam::hexRef8::hexRef8
|
|||
{
|
||||
FatalErrorIn
|
||||
(
|
||||
"hexRef8::hexRef8(const polyMesh&, const labelList&"
|
||||
", const labelList&, const refinementHistory&)"
|
||||
"hexRef8::hexRef8(const polyMesh&, const labelField&"
|
||||
", const labelField&, const refinementHistory&)"
|
||||
) << "Incorrect cellLevel or pointLevel size." << endl
|
||||
<< "Number of cells in mesh:" << mesh_.nCells()
|
||||
<< " does not equal size of cellLevel:" << cellLevel_.size() << endl
|
||||
|
@ -1876,8 +1876,8 @@ Foam::hexRef8::hexRef8
|
|||
Foam::hexRef8::hexRef8
|
||||
(
|
||||
const polyMesh& mesh,
|
||||
const labelList& cellLevel,
|
||||
const labelList& pointLevel
|
||||
const labelField& cellLevel,
|
||||
const labelField& pointLevel
|
||||
)
|
||||
:
|
||||
mesh_(mesh),
|
||||
|
|
|
@ -43,6 +43,7 @@ SourceFiles
|
|||
#include "removeFaces.H"
|
||||
#include "refinementHistory.H"
|
||||
#include "PackedList.H"
|
||||
#include "labelIOField.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
|
@ -68,10 +69,10 @@ class hexRef8
|
|||
const polyMesh& mesh_;
|
||||
|
||||
//- Per cell the refinement level
|
||||
labelIOList cellLevel_;
|
||||
labelIOField cellLevel_;
|
||||
|
||||
//- Per point the refinement level
|
||||
labelIOList pointLevel_;
|
||||
labelIOField pointLevel_;
|
||||
|
||||
//- Typical edge length between unrefined points
|
||||
const scalar level0Edge_;
|
||||
|
@ -326,8 +327,8 @@ public:
|
|||
hexRef8
|
||||
(
|
||||
const polyMesh& mesh,
|
||||
const labelList& cellLevel,
|
||||
const labelList& pointLevel,
|
||||
const labelField& cellLevel,
|
||||
const labelField& pointLevel,
|
||||
const refinementHistory& history
|
||||
);
|
||||
|
||||
|
@ -335,8 +336,8 @@ public:
|
|||
hexRef8
|
||||
(
|
||||
const polyMesh& mesh,
|
||||
const labelList& cellLevel,
|
||||
const labelList& pointLevel
|
||||
const labelField& cellLevel,
|
||||
const labelField& pointLevel
|
||||
);
|
||||
|
||||
|
||||
|
@ -344,12 +345,12 @@ public:
|
|||
|
||||
// Access
|
||||
|
||||
const labelIOList& cellLevel() const
|
||||
const labelIOField& cellLevel() const
|
||||
{
|
||||
return cellLevel_;
|
||||
}
|
||||
|
||||
const labelIOList& pointLevel() const
|
||||
const labelIOField& pointLevel() const
|
||||
{
|
||||
return pointLevel_;
|
||||
}
|
||||
|
|
|
@ -264,8 +264,8 @@ void Foam::multiDirRefinement::refineHex8
|
|||
hexRef8 hexRefiner
|
||||
(
|
||||
mesh,
|
||||
labelList(mesh.nCells(), 0), // cellLevel
|
||||
labelList(mesh.nPoints(), 0), // pointLevel
|
||||
labelField(mesh.nCells(), 0), // cellLevel
|
||||
labelField(mesh.nPoints(), 0), // pointLevel
|
||||
refinementHistory
|
||||
(
|
||||
IOobject
|
||||
|
|
Reference in a new issue