decompositionMethods::zoltan: Added support for processors contaning zero cells

This update allows Zoltan to be used by snappyHexMesh to redistribute the mesh
after refinement and sub-setting even if some processors loose all their cells
in the process.
This commit is contained in:
Henry Weller
2024-06-21 14:54:35 +01:00
parent 2f7185d73a
commit adfc80c412
2 changed files with 44 additions and 22 deletions

View File

@ -2,7 +2,7 @@
========= | ========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org \\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2023 OpenFOAM Foundation \\ / A nd | Copyright (C) 2023-2024 OpenFOAM Foundation
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License

View File

@ -26,6 +26,7 @@ License
#include "zoltan.H" #include "zoltan.H"
#include "globalIndex.H" #include "globalIndex.H"
#include "PstreamGlobals.H" #include "PstreamGlobals.H"
#include "Tuple3.H"
#include "addToRunTimeSelectionTable.H" #include "addToRunTimeSelectionTable.H"
#include "sigFpe.H" #include "sigFpe.H"
@ -81,28 +82,32 @@ static void get_vertex_list
int* ierr int* ierr
) )
{ {
const Foam::Tuple2<const Foam::pointField&, const Foam::scalarField&>& const Foam::Tuple3
vertexData = <
const Foam::pointField&,
const Foam::scalarField&,
const Foam::globalIndex&
>& vertexData =
*static_cast *static_cast
< <
const Foam::Tuple2 const Foam::Tuple3
< <
const Foam::pointField&, const Foam::pointField&,
const Foam::scalarField& const Foam::scalarField&,
const Foam::globalIndex&
>* >*
>(data); >(data);
const Foam::pointField& points = vertexData.first(); const Foam::pointField& points = vertexData.first();
const Foam::scalarField& weights = vertexData.second(); const Foam::scalarField& weights = vertexData.second();
const Foam::globalIndex& globalMap = vertexData.third();
if (wgt_dim != weights.size()/points.size()) if (points.size() && wgt_dim != weights.size()/points.size())
{ {
*ierr = ZOLTAN_FATAL; *ierr = ZOLTAN_FATAL;
return; return;
} }
Foam::globalIndex globalMap(points.size());
for (Foam::label i=0; i<points.size(); i++) for (Foam::label i=0; i<points.size(); i++)
{ {
localIDs[i] = i; localIDs[i] = i;
@ -216,13 +221,25 @@ static void get_edge_list
int* ierr int* ierr
) )
{ {
const Foam::CompactListList<Foam::label>& adjacency = const Foam::Tuple2
*static_cast<const Foam::CompactListList<Foam::label>*>(data); <
const Foam::CompactListList<Foam::label>&,
const Foam::globalIndex&
>& edgeData =
*static_cast
<
const Foam::Tuple2
<
const Foam::CompactListList<Foam::label>&,
const Foam::globalIndex&
>*
>(data);
const Foam::CompactListList<Foam::label>& adjacency = edgeData.first();
const Foam::globalIndex& globalMap = edgeData.second();
const Foam::labelList& m(adjacency.m()); const Foam::labelList& m(adjacency.m());
Foam::globalIndex globalMap(nPoints);
if if
( (
(nGID != 1) (nGID != 1)
@ -251,14 +268,6 @@ Foam::label Foam::decompositionMethods::zoltan::decompose
List<label>& decomp List<label>& decomp
) const ) const
{ {
if (Pstream::parRun() && !points.size())
{
Pout<< "No points on processor " << Pstream::myProcNo() << nl
<< " Zoltan cannot redistribute without points"
"on all processors." << endl;
FatalErrorInFunction << exit(FatalError);
}
stringList args(1); stringList args(1);
args[0] = "zoltan"; args[0] = "zoltan";
@ -330,10 +339,17 @@ Foam::label Foam::decompositionMethods::zoltan::decompose
<< exit(FatalIOError); << exit(FatalIOError);
} }
globalIndex globalMap(points.size());
void* pointsPtr = &const_cast<pointField&>(points); void* pointsPtr = &const_cast<pointField&>(points);
Zoltan_Set_Num_Obj_Fn(zz, get_number_of_vertices, pointsPtr); Zoltan_Set_Num_Obj_Fn(zz, get_number_of_vertices, pointsPtr);
Tuple2<const pointField&, const scalarField&> vertexData(points, pWeights); Tuple3<const pointField&, const scalarField&, const globalIndex&> vertexData
(
points,
pWeights,
globalMap
);
Zoltan_Set_Obj_List_Fn(zz, get_vertex_list, &vertexData); Zoltan_Set_Obj_List_Fn(zz, get_vertex_list, &vertexData);
// Callbacks for geometry // Callbacks for geometry
@ -343,7 +359,13 @@ Foam::label Foam::decompositionMethods::zoltan::decompose
// Callbacks for connectivity // Callbacks for connectivity
void* adjacencyPtr = &const_cast<CompactListList<label>&>(adjacency); void* adjacencyPtr = &const_cast<CompactListList<label>&>(adjacency);
Zoltan_Set_Num_Edges_Multi_Fn(zz, get_num_edges_list, adjacencyPtr); Zoltan_Set_Num_Edges_Multi_Fn(zz, get_num_edges_list, adjacencyPtr);
Zoltan_Set_Edge_List_Multi_Fn(zz, get_edge_list, adjacencyPtr);
Tuple2<const CompactListList<label>&, const globalIndex&> edgeData
(
adjacency,
globalMap
);
Zoltan_Set_Edge_List_Multi_Fn(zz, get_edge_list, &edgeData);
int changed, num_imp, num_exp, *imp_procs, *exp_procs; int changed, num_imp, num_exp, *imp_procs, *exp_procs;
int *imp_to_part, *exp_to_part; int *imp_to_part, *exp_to_part;