45#ifndef _ZOLTAN2_ALGZOLTAN_HPP_
46#define _ZOLTAN2_ALGZOLTAN_HPP_
57#include <zoltan_cpp.h>
58#include <zoltan_partition_tree.h>
76template <
typename Adapter>
84 typedef typename Adapter::scalar_t scalar_t;
87 typedef typename Adapter::userCoord_t userCoord_t;
89 const RCP<const Environment> env;
90 const RCP<const Comm<int> > problemComm;
91 const RCP<const typename Adapter::base_adapter_t> adapter;
92 RCP<const Model<Adapter> > model;
97 void setMPIComm(
const RCP<
const Comm<int> > &problemComm__) {
98# ifdef HAVE_ZOLTAN2_MPI
99 mpicomm = Teuchos::getRawMpiComm(*problemComm__);
101 mpicomm = MPI_COMM_WORLD;
110 Zoltan_Initialize(argc, argv, &ver);
113 void setCallbacksIDs()
115 zz->Set_Num_Obj_Fn(zoltanNumObj<Adapter>, (
void *) &(*adapter));
116 zz->Set_Obj_List_Fn(zoltanObjList<Adapter>, (
void *) &(*adapter));
118 const part_t *myparts;
119 adapter->getPartsView(myparts);
121 zz->Set_Part_Multi_Fn(zoltanParts<Adapter>, (
void *) &(*adapter));
125 zz->Set_Param(
"NUM_GID_ENTRIES", tmp);
127 zz->Set_Param(
"NUM_LID_ENTRIES", tmp);
130 template <
typename AdapterWithCoords>
131 void setCallbacksGeom(
const AdapterWithCoords *ia)
136 zz->Set_Num_Geom_Fn(zoltanNumGeom<AdapterWithCoords>, (
void *) ia);
137 zz->Set_Geom_Multi_Fn(zoltanGeom<AdapterWithCoords>, (
void *) ia);
140 void setCallbacksGraph(
147 void setCallbacksGraph(
154 void setCallbacksGraph(
161 void setCallbacksHypergraph(
167 zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withMatrixAdapter<Adapter>,
169 zz->Set_HG_CS_Fn(zoltanHGCS_withMatrixAdapter<Adapter>,
178 void setCallbacksHypergraph(
181 zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withGraphAdapter<Adapter>,
183 zz->Set_HG_CS_Fn(zoltanHGCS_withGraphAdapter<Adapter>,
186 if (adp->getNumWeightsPerEdge() != 0) {
187 if (adp->getNumWeightsPerEdge() > 1) {
188 std::cout <<
"Zoltan2 warning: getNumWeightsPerEdge() returned "
189 << adp->getNumWeightsPerEdge() <<
" but PHG supports only "
190 <<
" one weight per edge; only first weight will be used."
193 zz->Set_HG_Size_Edge_Wts_Fn(zoltanHGSizeEdgeWts_withGraphAdapter<Adapter>,
194 (
void *) &(*adapter));
195 zz->Set_HG_Edge_Wts_Fn(zoltanHGEdgeWts_withGraphAdapter<Adapter>,
196 (
void *) &(*adapter));
203 const Teuchos::ParameterList &pl = env->getParameters();
205 const Teuchos::ParameterEntry *pe = pl.getEntryPtr(
"hypergraph_model_type");
206 std::string model_type(
"traditional");
208 model_type = pe->getValue<std::string>(&model_type);
211 if (model_type==
"ghosting" ||
212 !adp->areEntityIDsUnique(adp->getPrimaryEntityType())) {
219 zz->Set_Num_Obj_Fn(zoltanHGNumObj_withModel<Adapter>, (
void *) &(*mdl));
220 zz->Set_Obj_List_Fn(zoltanHGObjList_withModel<Adapter>, (
void *) &(*mdl));
222 zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withModel<Adapter>, (
void *) &(*mdl));
223 zz->Set_HG_CS_Fn(zoltanHGCS_withModel<Adapter>, (
void *) &(*mdl));
227 zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withMeshAdapter<Adapter>,
229 zz->Set_HG_CS_Fn(zoltanHGCS_withMeshAdapter<Adapter>,
239 virtual bool isPartitioningTreeBinary()
const
245 void rcb_recursive_partitionTree_calculations(
248 std::vector<part_t> & splitRangeBeg,
249 std::vector<part_t> & splitRangeEnd)
const
258 int err = Zoltan_RCB_Partition_Tree(zz->Get_C_Handle(),
259 arrayIndex - numParts + 1,
260 &parent, &left_leaf, &right_leaf);
262 throw std::logic_error(
"Zoltan_RCB_Partition_Tree returned in error." );
267 rcb_recursive_partitionTree_calculations(left_leaf+numParts-1, numParts,
268 splitRangeBeg, splitRangeEnd);
271 rcb_recursive_partitionTree_calculations(right_leaf+numParts-1, numParts,
272 splitRangeBeg, splitRangeEnd);
277 int leftIndex = (left_leaf > 0) ? (left_leaf-1+numParts) : (-left_leaf);
278 int rightIndex = (right_leaf > 0) ? (right_leaf-1+numParts) : (-right_leaf);
279 splitRangeBeg[arrayIndex] = splitRangeBeg[leftIndex];
280 splitRangeEnd[arrayIndex] = splitRangeEnd[rightIndex];
283 if(splitRangeBeg[rightIndex] != splitRangeEnd[leftIndex]) {
284 throw std::logic_error(
"RCB expected left_leaf indices and right leaf"
285 " indices to be continuous but it was not so." );
290 void rcb_getPartitionTree(part_t numParts,
291 part_t & numTreeVerts,
292 std::vector<part_t> & permPartNums,
293 std::vector<part_t> & splitRangeBeg,
294 std::vector<part_t> & splitRangeEnd,
295 std::vector<part_t> & treeVertParents)
const
303 part_t numTreeNodes = 2 * numParts - 1;
304 numTreeVerts = numTreeNodes - 1;
306 permPartNums.resize(numParts);
307 for(part_t n = 0; n < numParts; ++n) {
311 treeVertParents.resize(numTreeNodes);
319 part_t saveRootNodeChildrenA = -1;
320 part_t saveRootNodeChildrenB = -1;
321 part_t saveFinalNodeChildrenA = -1;
322 part_t saveFinalNodeChildrenB = -1;
323 for(part_t n = numParts; n < numTreeNodes; ++n) {
327 int err = Zoltan_RCB_Partition_Tree(zz->Get_C_Handle(),
328 static_cast<int>(n - numParts) + 1,
329 &parent, &left_leaf, &right_leaf);
331 throw std::logic_error(
"Zoltan_RCB_Partition_Tree returned in error.");
333 part_t leftIndex = (left_leaf > 0) ? (left_leaf-1+numParts) : (-left_leaf);
334 part_t rightIndex = (right_leaf > 0) ? (right_leaf-1+numParts) : (-right_leaf);
335 treeVertParents[leftIndex] = n;
336 treeVertParents[rightIndex] = n;
338 if(parent == 1 || parent == -1) {
340 saveRootNodeChildrenA = leftIndex;
341 saveRootNodeChildrenB = rightIndex;
343 if(n == numTreeNodes-1) {
344 saveFinalNodeChildrenA = leftIndex;
345 saveFinalNodeChildrenB = rightIndex;
348 treeVertParents[rootNode] = -1;
350 splitRangeBeg.resize(numTreeNodes);
351 splitRangeEnd.resize(numTreeNodes);
353 for(part_t n = 0; n < numParts; ++n) {
354 splitRangeBeg[n] = n;
355 splitRangeEnd[n] = n + 1;
361 rcb_recursive_partitionTree_calculations(rootNode, numParts, splitRangeBeg,
365 std::swap(treeVertParents[rootNode], treeVertParents[numTreeNodes-1]);
367 treeVertParents[saveFinalNodeChildrenA] = rootNode;
368 treeVertParents[saveFinalNodeChildrenB] = rootNode;
370 if(saveRootNodeChildrenA == numTreeNodes - 1) {
371 saveRootNodeChildrenA = rootNode;
373 if(saveRootNodeChildrenB == numTreeNodes - 1) {
374 saveRootNodeChildrenB = rootNode;
376 treeVertParents[saveRootNodeChildrenA] = numTreeNodes - 1;
377 treeVertParents[saveRootNodeChildrenB] = numTreeNodes - 1;
379 std::swap(splitRangeBeg[rootNode], splitRangeBeg[numTreeNodes-1]);
380 std::swap(splitRangeEnd[rootNode], splitRangeEnd[numTreeNodes-1]);
385 void phg_getPartitionTree(part_t numParts,
386 part_t & numTreeVerts,
387 std::vector<part_t> & permPartNums,
388 std::vector<part_t> & splitRangeBeg,
389 std::vector<part_t> & splitRangeEnd,
390 std::vector<part_t> & treeVertParents)
const
398 int tree_array_size = -1;
399 int err = Zoltan_PHG_Partition_Tree_Size(
400 zz->Get_C_Handle(), &tree_array_size);
402 throw std::logic_error(
"Zoltan_PHG_Partition_Tree_Size returned error.");
413 part_t numTreeNodes = 0;
414 std::vector<part_t> mapIndex(tree_array_size);
415 part_t trackNonTerminalIndex = numParts;
416 for(part_t n = 0; n < static_cast<part_t>(tree_array_size); ++n) {
417 part_t phgIndex = n + 1;
420 err = Zoltan_PHG_Partition_Tree(
421 zz->Get_C_Handle(), phgIndex, &lo_index, &hi_index);
425 mapIndex[n] = (lo_index == hi_index) ?
427 (trackNonTerminalIndex++);
432 mapIndex[0] = numTreeNodes - 1;
434 numTreeVerts = numTreeNodes - 1;
436 permPartNums.resize(numParts);
437 for(part_t n = 0; n < numParts; ++n) {
443 treeVertParents.resize(numTreeNodes);
444 splitRangeBeg.resize(numTreeNodes);
445 splitRangeEnd.resize(numTreeNodes);
447 for(part_t n = 0; n < tree_array_size; ++n) {
448 part_t phgIndex = n + 1;
451 err = Zoltan_PHG_Partition_Tree(
452 zz->Get_C_Handle(), phgIndex, &lo_index, &hi_index);
454 throw std::logic_error(
"Zoltan_PHG_Partition_Tree returned in error.");
458 part_t finalIndex = mapIndex[n];
467 part_t parentPHGIndex =
468 ((phgIndex%2) == 0) ? (phgIndex/2) : ((phgIndex-1)/2);
471 treeVertParents[finalIndex] = (n==0) ? -1 : mapIndex[parentPHGIndex-1];
473 splitRangeBeg[finalIndex] =
static_cast<part_t
>(lo_index);
474 splitRangeEnd[finalIndex] =
static_cast<part_t
>(hi_index+1);
480 void getPartitionTree(part_t numParts,
481 part_t & numTreeVerts,
482 std::vector<part_t> & permPartNums,
483 std::vector<part_t> & splitRangeBeg,
484 std::vector<part_t> & splitRangeEnd,
485 std::vector<part_t> & treeVertParents)
const
488 const Teuchos::ParameterList &pl = env->getParameters();
489 bool keep_partition_tree =
false;
490 const Teuchos::ParameterEntry * pe = pl.getEntryPtr(
"keep_partition_tree");
492 keep_partition_tree = pe->getValue(&keep_partition_tree);
493 if(!keep_partition_tree) {
494 throw std::logic_error(
495 "Requested tree when param keep_partition_tree is off.");
501 const Teuchos::ParameterList & zoltan_pl = pl.sublist(
"zoltan_parameters");
502 std::string lb_method;
503 pe = zoltan_pl.getEntryPtr(
"LB_METHOD");
505 lb_method = pe->getValue(&lb_method);
507 if(lb_method ==
"phg") {
508 phg_getPartitionTree(numParts, numTreeVerts, permPartNums,
509 splitRangeBeg, splitRangeEnd, treeVertParents);
511 else if(lb_method ==
"rcb") {
512 rcb_getPartitionTree(numParts, numTreeVerts, permPartNums,
513 splitRangeBeg, splitRangeEnd, treeVertParents);
516 throw std::logic_error(
"Did not recognize LB_METHOD: " + lb_method);
528 const RCP<
const Comm<int> > &problemComm__,
530 env(env__), problemComm(problemComm__), adapter(adapter__)
532 setMPIComm(problemComm__);
534 zz = rcp(
new Zoltan(mpicomm));
539 const RCP<
const Comm<int> > &problemComm__,
541 env(env__), problemComm(problemComm__), adapter(adapter__)
543 setMPIComm(problemComm__);
545 zz = rcp(
new Zoltan(mpicomm));
547 setCallbacksGeom(&(*adapter));
551 const RCP<
const Comm<int> > &problemComm__,
553 env(env__), problemComm(problemComm__), adapter(adapter__)
555 setMPIComm(problemComm__);
557 zz = rcp(
new Zoltan(mpicomm));
559 setCallbacksGraph(adapter);
560 setCallbacksHypergraph(adapter);
561 if (adapter->coordinatesAvailable()) {
562 setCallbacksGeom(adapter->getCoordinateInput());
567 const RCP<
const Comm<int> > &problemComm__,
569 env(env__), problemComm(problemComm__), adapter(adapter__)
571 setMPIComm(problemComm__);
573 zz = rcp(
new Zoltan(mpicomm));
575 setCallbacksGraph(adapter);
576 setCallbacksHypergraph(adapter);
577 if (adapter->coordinatesAvailable()) {
578 setCallbacksGeom(adapter->getCoordinateInput());
583 const RCP<
const Comm<int> > &problemComm__,
585 env(env__), problemComm(problemComm__), adapter(adapter__)
587 setMPIComm(problemComm__);
589 zz = rcp(
new Zoltan(mpicomm));
591 setCallbacksGraph(adapter);
595 setCallbacksHypergraph(adapter);
596 setCallbacksGeom(&(*adapter));
604template <
typename Adapter>
612 size_t numGlobalParts = solution->getTargetGlobalNumberOfParts();
614 sprintf(paramstr,
"%lu", numGlobalParts);
615 zz->Set_Param(
"NUM_GLOBAL_PARTS", paramstr);
617 int wdim = adapter->getNumWeightsPerID();
618 sprintf(paramstr,
"%d", wdim);
619 zz->Set_Param(
"OBJ_WEIGHT_DIM", paramstr);
621 const Teuchos::ParameterList &pl = env->getParameters();
624 const Teuchos::ParameterEntry *pe = pl.getEntryPtr(
"imbalance_tolerance");
627 tolerance = pe->getValue<
double>(&tolerance);
628 sprintf(str,
"%f", tolerance);
629 zz->Set_Param(
"IMBALANCE_TOL", str);
632 pe = pl.getEntryPtr(
"partitioning_approach");
634 std::string approach;
635 approach = pe->getValue<std::string>(&approach);
636 if (approach ==
"partition")
637 zz->Set_Param(
"LB_APPROACH",
"PARTITION");
639 zz->Set_Param(
"LB_APPROACH",
"REPARTITION");
642 pe = pl.getEntryPtr(
"partitioning_objective");
644 std::string strChoice = pe->getValue<std::string>(&strChoice);
645 if (strChoice == std::string(
"multicriteria_minimize_total_weight"))
646 zz->Set_Param(
"RCB_MULTICRITERIA_NORM",
"1");
647 else if (strChoice == std::string(
"multicriteria_balance_total_maximum"))
648 zz->Set_Param(
"RCB_MULTICRITERIA_NORM",
"2");
649 else if (strChoice == std::string(
"multicriteria_minimize_maximum_weight"))
650 zz->Set_Param(
"RCB_MULTICRITERIA_NORM",
"3");
657 bool keep_partition_tree =
false;
658 pe = pl.getEntryPtr(
"keep_partition_tree");
660 keep_partition_tree = pe->getValue(&keep_partition_tree);
661 if (keep_partition_tree) {
665 zz->Set_Param(
"KEEP_CUTS",
"1");
666 zz->Set_Param(
"PHG_KEEP_TREE",
"1");
670 pe = pl.getEntryPtr(
"rectilinear");
672 bool val = pe->getValue(&val);
674 zz->Set_Param(
"RCB_RECTILINEAR_BLOCKS",
"1");
679 const Teuchos::ParameterList &zpl = pl.sublist(
"zoltan_parameters");
680 for (ParameterList::ConstIterator iter = zpl.begin();
681 iter != zpl.end(); iter++) {
682 const std::string &zname = pl.name(iter);
684 std::string zval = pl.entry(iter).getValue(&zval);
685 zz->Set_Param(zname.c_str(), zval.c_str());
688 catch (std::exception &e) {
693 int pdim = (wdim > 1 ? wdim : 1);
694 for (
int d = 0; d < pdim; d++) {
695 if (!solution->criteriaHasUniformPartSizes(d)) {
696 float *partsizes =
new float[numGlobalParts];
697 int *partidx =
new int[numGlobalParts];
698 int *wgtidx =
new int[numGlobalParts];
699 for (
size_t i=0; i<numGlobalParts; i++) partidx[i] = i;
700 for (
size_t i=0; i<numGlobalParts; i++) wgtidx[i] = d;
701 for (
size_t i=0; i<numGlobalParts; i++)
702 partsizes[i] = solution->getCriteriaPartSize(0, i);
703 zz->LB_Set_Part_Sizes(1, numGlobalParts, partidx, wgtidx, partsizes);
716 ZOLTAN_ID_PTR dGids = NULL, dLids = NULL;
717 int *dProcs = NULL, *dParts = NULL;
719 ZOLTAN_ID_PTR oGids = NULL, oLids = NULL;
720 int *oProcs = NULL, *oParts = NULL;
722 zz->Set_Param(
"RETURN_LISTS",
"PARTS");
725 int ierr = zz->LB_Partition(changed, nGidEnt, nLidEnt,
726 nDummy, dGids, dLids, dProcs, dParts,
727 nObj, oGids, oLids, oProcs, oParts);
729 env->globalInputAssertion(__FILE__, __LINE__,
"Zoltan LB_Partition",
738 numObjects=model->getLocalNumObjects();
742 ArrayRCP<part_t> partList(
new part_t[numObjects], 0, numObjects,
true);
743 for (
int i = 0; i < nObj; i++) {
746 partList[tmp] = oParts[i];
758 Teuchos::RCP<const map_t> mapWithCopies;
759 Teuchos::RCP<const map_t> oneToOneMap;
762 typedef Tpetra::Vector<scalar_t, lno_t, gno_t> vector_t;
763 vector_t vecWithCopies(mapWithCopies);
764 vector_t oneToOneVec(oneToOneMap);
767 assert(nObj ==
lno_t(oneToOneMap->getNodeNumElements()));
768 for (lno_t i = 0; i < nObj; i++)
769 oneToOneVec.replaceLocalValue(i, oParts[i]);
772 Teuchos::RCP<const Tpetra::Import<lno_t, gno_t> > importer =
773 Tpetra::createImport<lno_t, gno_t>(oneToOneMap, mapWithCopies);
774 vecWithCopies.doImport(oneToOneVec, *importer, Tpetra::REPLACE);
777 lno_t nlocal =
lno_t(mapWithCopies->getNodeNumElements());
778 for (lno_t i = 0; i < nlocal; i++)
779 partList[i] = vecWithCopies.getData()[i];
782 solution->setParts(partList);
785 zz->LB_Free_Part(&oGids, &oLids, &oProcs, &oParts);
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
callback functions for the Zoltan package (templated on Adapter)
Defines the Model interface.
Defines the PartitioningSolution class.
Gathering definitions used in software development.
Traits class to handle conversions between gno_t/lno_t and TPL data types (e.g., ParMETIS's idx_t,...
A gathering of useful namespace methods.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const MeshAdapter< user_t > > &adapter__)
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const IdentifierAdapter< user_t > > &adapter__)
void partition(const RCP< PartitioningSolution< Adapter > > &solution)
Partitioning method.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const MatrixAdapter< user_t, userCoord_t > > &adapter__)
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const VectorAdapter< user_t > > &adapter__)
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const GraphAdapter< user_t, userCoord_t > > &adapter__)
Algorithm defines the base class for all algorithms.
GraphAdapter defines the interface for graph-based user data.
HyperGraphModel defines the interface required for hyper graph models.
void getVertexMaps(Teuchos::RCP< const map_t > &copiesMap, Teuchos::RCP< const map_t > &onetooneMap) const
Sets pointers to the vertex map with copies and the vertex map without copies Note: the pointers will...
IdentifierAdapter defines the interface for identifiers.
MatrixAdapter defines the adapter interface for matrices.
MeshAdapter defines the interface for mesh input.
The base class for all model classes.
A PartitioningSolution is a solution to a partitioning problem.
VectorAdapter defines the interface for vector input.
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
Created by mbenlioglu on Aug 31, 2020.
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
@ BASIC_ASSERTION
fast typical checks for valid arguments
SparseMatrixAdapter_t::part_t part_t
static void ASSIGN(first_t &a, second_t b)