46 namespace decompositionMethods
78 ZOLTAN_ID_PTR globalIDs,
79 ZOLTAN_ID_PTR localIDs,
107 *ierr = ZOLTAN_FATAL;
114 globalIDs[i] = globalMap.
toGlobal(i);
116 for(
int j=0; j<wgt_dim; j++)
118 obj_wgts[wgt_dim*i + j] = weights[wgt_dim*i + j];
138 ZOLTAN_ID_PTR globalIDs,
139 ZOLTAN_ID_PTR localIDs,
156 *ierr = ZOLTAN_FATAL;
182 ZOLTAN_ID_PTR globalIDs,
183 ZOLTAN_ID_PTR localIDs,
193 if ((nGID != 1) || (nLID != 1) || (
nPoints != sizes.
size()))
195 *ierr = ZOLTAN_FATAL;
201 numEdges[i] = sizes[i];
214 ZOLTAN_ID_PTR globalIDs,
215 ZOLTAN_ID_PTR localIDs,
217 ZOLTAN_ID_PTR nborGID,
249 *ierr = ZOLTAN_FATAL;
265 const CompactListList<label>& adjacency,
276 for (
label i = 0; i < argc; i++)
278 argv[i] = strdup(
args[i].c_str());
282 int rc = Zoltan_Initialize(argc, argv, &ver);
295 Zoltan_Set_Param(zz,
"return_lists",
"export");
296 Zoltan_Set_Param(zz,
"obj_weight_dim",
name(
nWeights).c_str());
297 Zoltan_Set_Param(zz,
"edge_weight_dim",
"0");
300 Zoltan_Set_Param(zz,
"debug_level",
"0");
301 Zoltan_Set_Param(zz,
"imbalance_tol",
"1.05");
303 word lb_method(
"graph");
306 Zoltan_Set_Param(zz,
"lb_method", lb_method.c_str());
307 Zoltan_Set_Param(zz,
"lb_approach",
"repartition");
311 const dictionary& coeffsDict_ =
318 if (!iter().isDict())
320 const word& key = iter().keyword();
321 const word value(iter().stream());
325 Info<< typeName <<
" : setting parameter " << key
326 <<
" to " << value <<
endl;
329 Zoltan_Set_Param(zz, key.c_str(), value.c_str());
334 if (
nWeights > 1 && lb_method !=
"rcb")
337 <<
"Multiple constraints specified for lb_method " << lb_method
338 <<
" but is only supported by the rcb method"
342 globalIndex globalMap(
points.size());
347 Tuple3<const pointField&, const scalarField&, const globalIndex&> vertexData
360 void* adjacencyPtr = &
const_cast<CompactListList<label>&
>(adjacency);
363 Tuple2<const CompactListList<label>&,
const globalIndex&> edgeData
370 int changed, num_imp, num_exp, *imp_procs, *exp_procs;
371 int *imp_to_part, *exp_to_part;
372 int num_gid_entries, num_lid_entries;
373 ZOLTAN_ID_PTR imp_global_ids, exp_global_ids;
374 ZOLTAN_ID_PTR imp_local_ids, exp_local_ids;
379 int oldExcepts = fedisableexcept
408 feenableexcept(oldExcepts);
411 if (debug && changed)
413 Pout <<
"Zoltan number to move " << changed <<
" " << num_exp <<
endl;
416 decomp.setSize(
points.size());
421 for(
int i=0; i<num_exp; i++)
423 decomp[exp_local_ids[i]] = exp_procs[i];
431 nProcCells[decomp[i]]++;
434 reduce(nProcCells, ListOp<sumOp<label>>());
440 Pout<<
" No cells allocated to this processor"
441 ", keeping first cell"
468 for (
label i = 0; i < argc; i++)
497 <<
"Can use this decomposition method only for the whole mesh"
499 <<
"and supply one coordinate (cellCentre) for every cell." <<
endl
500 <<
"The number of coordinates " << cellCentres.
size() <<
endl
501 <<
"The number of cells in the mesh " << mesh.
nCells()
542 <<
"Size of cell-to-coarse map " << agglom.
size()
543 <<
" differs from number of cells in mesh " << mesh.
nCells()
571 forAll(fineDistribution, i)
573 fineDistribution[i] = decomp[agglom[i]];
576 return fineDistribution;
587 if (
points.size() != globalPointPoints.
size())
590 <<
"Inconsistent number of cells (" << globalPointPoints.
size()
591 <<
") and number of cell centres (" <<
points.size()
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Macros for easy insertion into run-time selection tables.
A packed storage unstructured matrix of objects of type <T> using an offset table for access.
void size(const label)
Override size to be inconsistent with allocated storage.
A 2-tuple for storing two objects of different types.
A 3-tuple for storing three objects of different types.
label size() const
Return the primary size, i.e. the number of rows.
labelList sizes() const
Return sizes (to be used e.g. for construction)
static label nProcs(const label communicator=0)
Number of processes in parallel run.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
static const direction nComponents
Number of components in this vector space.
label size() const
Return the number of arguments.
Abstract base class for decomposition.
label nWeights(const pointField &points, const scalarField &pointWeights) const
Return the number of weights per point.
dictionary decompositionDict_
Zoltan redistribution in parallel.
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Inherit decompose from decompositionMethod.
zoltan(const dictionary &decompositionDict)
Construct given the decomposition dictionary and mesh.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
label whichProcID(const label i) const
Which processor does global come from? Binary search.
label toGlobal(const label i) const
From local to global.
Mesh consisting of general polyhedral cells.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
volScalarField scalarField(fieldObject, mesh)
defineTypeNameAndDebug(metis, 0)
addToRunTimeSelectionTable(decompositionMethod, metis, decomposer)
errorManipArg< error, int > exit(error &err, const int errNo=1)
pointField vertices(const blockVertexList &bvl)
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
vectorField pointField
pointField is a vectorField.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
prefixOSstream Pout(cout, "Pout")
List< string > stringList
A List of strings.
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Foam::argList args(argc, argv)
static int get_number_of_vertices(void *data, int *ierr)
static void get_num_edges_list(void *data, int nGID, int nLID, int nPoints, ZOLTAN_ID_PTR globalIDs, ZOLTAN_ID_PTR localIDs, int *numEdges, int *ierr)
static void get_geom_list(void *data, int nGID, int nLID, int nPoints, ZOLTAN_ID_PTR globalIDs, ZOLTAN_ID_PTR localIDs, int nDim, double *vertices, int *ierr)
static int get_mesh_dim(void *data, int *ierr)
static void get_vertex_list(void *data, int nGID, int nLID, ZOLTAN_ID_PTR globalIDs, ZOLTAN_ID_PTR localIDs, int wgt_dim, float *obj_wgts, int *ierr)
static void get_edge_list(void *data, int nGID, int nLID, int nPoints, ZOLTAN_ID_PTR globalIDs, ZOLTAN_ID_PTR localIDs, int *num_edges, ZOLTAN_ID_PTR nborGID, int *nborProc, int wgt_dim, float *ewgts, int *ierr)