45 #ifndef _ZOLTAN2_ALGZOLTAN_HPP_ 46 #define _ZOLTAN2_ALGZOLTAN_HPP_ 57 #include <zoltan_cpp.h> 58 #include <zoltan_partition_tree.h> 76 template <
typename Adapter>
82 typedef typename Adapter::lno_t lno_t;
83 typedef typename Adapter::gno_t gno_t;
84 typedef typename Adapter::scalar_t scalar_t;
86 typedef typename Adapter::user_t user_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>,
181 const Teuchos::ParameterList &pl = env->getParameters();
183 const Teuchos::ParameterEntry *pe = pl.getEntryPtr(
"hypergraph_model_type");
184 std::string model_type(
"traditional");
186 model_type = pe->getValue<std::string>(&model_type);
189 if (model_type==
"ghosting" ||
190 !adp->areEntityIDsUnique(adp->getPrimaryEntityType())) {
197 zz->Set_Num_Obj_Fn(zoltanHGNumObj_withModel<Adapter>, (
void *) &(*mdl));
198 zz->Set_Obj_List_Fn(zoltanHGObjList_withModel<Adapter>, (
void *) &(*mdl));
200 zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withModel<Adapter>, (
void *) &(*mdl));
201 zz->Set_HG_CS_Fn(zoltanHGCS_withModel<Adapter>, (
void *) &(*mdl));
205 zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withMeshAdapter<Adapter>,
207 zz->Set_HG_CS_Fn(zoltanHGCS_withMeshAdapter<Adapter>,
217 virtual bool isPartitioningTreeBinary()
const 223 void rcb_recursive_partitionTree_calculations(
226 std::vector<part_t> & splitRangeBeg,
227 std::vector<part_t> & splitRangeEnd)
const 236 int err = Zoltan_RCB_Partition_Tree(zz->Get_C_Handle(),
237 arrayIndex - numParts + 1,
238 &parent, &left_leaf, &right_leaf);
240 throw std::logic_error(
"Zoltan_RCB_Partition_Tree returned in error." );
245 rcb_recursive_partitionTree_calculations(left_leaf+numParts-1, numParts,
246 splitRangeBeg, splitRangeEnd);
249 rcb_recursive_partitionTree_calculations(right_leaf+numParts-1, numParts,
250 splitRangeBeg, splitRangeEnd);
255 int leftIndex = (left_leaf > 0) ? (left_leaf-1+numParts) : (-left_leaf);
256 int rightIndex = (right_leaf > 0) ? (right_leaf-1+numParts) : (-right_leaf);
257 splitRangeBeg[arrayIndex] = splitRangeBeg[leftIndex];
258 splitRangeEnd[arrayIndex] = splitRangeEnd[rightIndex];
261 if(splitRangeBeg[rightIndex] != splitRangeEnd[leftIndex]) {
262 throw std::logic_error(
"RCB expected left_leaf indices and right leaf" 263 " indices to be continuous but it was not so." );
268 void rcb_getPartitionTree(part_t numParts,
269 part_t & numTreeVerts,
270 std::vector<part_t> & permPartNums,
271 std::vector<part_t> & splitRangeBeg,
272 std::vector<part_t> & splitRangeEnd,
273 std::vector<part_t> & treeVertParents)
const 281 part_t numTreeNodes = 2 * numParts - 1;
282 numTreeVerts = numTreeNodes - 1;
284 permPartNums.resize(numParts);
285 for(part_t n = 0; n < numParts; ++n) {
289 treeVertParents.resize(numTreeNodes);
297 part_t saveRootNodeChildrenA = -1;
298 part_t saveRootNodeChildrenB = -1;
299 part_t saveFinalNodeChildrenA = -1;
300 part_t saveFinalNodeChildrenB = -1;
301 for(part_t n = numParts; n < numTreeNodes; ++n) {
305 int err = Zoltan_RCB_Partition_Tree(zz->Get_C_Handle(),
306 static_cast<int>(n - numParts) + 1,
307 &parent, &left_leaf, &right_leaf);
309 throw std::logic_error(
"Zoltan_RCB_Partition_Tree returned in error.");
311 part_t leftIndex = (left_leaf > 0) ? (left_leaf-1+numParts) : (-left_leaf);
312 part_t rightIndex = (right_leaf > 0) ? (right_leaf-1+numParts) : (-right_leaf);
313 treeVertParents[leftIndex] = n;
314 treeVertParents[rightIndex] = n;
316 if(parent == 1 || parent == -1) {
318 saveRootNodeChildrenA = leftIndex;
319 saveRootNodeChildrenB = rightIndex;
321 if(n == numTreeNodes-1) {
322 saveFinalNodeChildrenA = leftIndex;
323 saveFinalNodeChildrenB = rightIndex;
326 treeVertParents[rootNode] = -1;
328 splitRangeBeg.resize(numTreeNodes);
329 splitRangeEnd.resize(numTreeNodes);
331 for(part_t n = 0; n < numParts; ++n) {
332 splitRangeBeg[n] = n;
333 splitRangeEnd[n] = n + 1;
339 rcb_recursive_partitionTree_calculations(rootNode, numParts, splitRangeBeg,
343 std::swap(treeVertParents[rootNode], treeVertParents[numTreeNodes-1]);
345 treeVertParents[saveFinalNodeChildrenA] = rootNode;
346 treeVertParents[saveFinalNodeChildrenB] = rootNode;
348 if(saveRootNodeChildrenA == numTreeNodes - 1) {
349 saveRootNodeChildrenA = rootNode;
351 if(saveRootNodeChildrenB == numTreeNodes - 1) {
352 saveRootNodeChildrenB = rootNode;
354 treeVertParents[saveRootNodeChildrenA] = numTreeNodes - 1;
355 treeVertParents[saveRootNodeChildrenB] = numTreeNodes - 1;
357 std::swap(splitRangeBeg[rootNode], splitRangeBeg[numTreeNodes-1]);
358 std::swap(splitRangeEnd[rootNode], splitRangeEnd[numTreeNodes-1]);
363 void phg_getPartitionTree(part_t numParts,
364 part_t & numTreeVerts,
365 std::vector<part_t> & permPartNums,
366 std::vector<part_t> & splitRangeBeg,
367 std::vector<part_t> & splitRangeEnd,
368 std::vector<part_t> & treeVertParents)
const 376 int tree_array_size = -1;
377 int err = Zoltan_PHG_Partition_Tree_Size(
378 zz->Get_C_Handle(), &tree_array_size);
380 throw std::logic_error(
"Zoltan_PHG_Partition_Tree_Size returned error.");
391 part_t numTreeNodes = 0;
392 std::vector<part_t> mapIndex(tree_array_size);
393 part_t trackNonTerminalIndex = numParts;
394 for(part_t n = 0; n < static_cast<part_t>(tree_array_size); ++n) {
395 part_t phgIndex = n + 1;
398 err = Zoltan_PHG_Partition_Tree(
399 zz->Get_C_Handle(), phgIndex, &lo_index, &hi_index);
403 mapIndex[n] = (lo_index == hi_index) ?
405 (trackNonTerminalIndex++);
410 mapIndex[0] = numTreeNodes - 1;
412 numTreeVerts = numTreeNodes - 1;
414 permPartNums.resize(numParts);
415 for(part_t n = 0; n < numParts; ++n) {
421 treeVertParents.resize(numTreeNodes);
422 splitRangeBeg.resize(numTreeNodes);
423 splitRangeEnd.resize(numTreeNodes);
425 for(part_t n = 0; n < tree_array_size; ++n) {
426 part_t phgIndex = n + 1;
429 err = Zoltan_PHG_Partition_Tree(
430 zz->Get_C_Handle(), phgIndex, &lo_index, &hi_index);
432 throw std::logic_error(
"Zoltan_PHG_Partition_Tree returned in error.");
436 part_t finalIndex = mapIndex[n];
445 part_t parentPHGIndex =
446 ((phgIndex%2) == 0) ? (phgIndex/2) : ((phgIndex-1)/2);
449 treeVertParents[finalIndex] = (n==0) ? -1 : mapIndex[parentPHGIndex-1];
451 splitRangeBeg[finalIndex] =
static_cast<part_t
>(lo_index);
452 splitRangeEnd[finalIndex] =
static_cast<part_t
>(hi_index+1);
458 void getPartitionTree(part_t numParts,
459 part_t & numTreeVerts,
460 std::vector<part_t> & permPartNums,
461 std::vector<part_t> & splitRangeBeg,
462 std::vector<part_t> & splitRangeEnd,
463 std::vector<part_t> & treeVertParents)
const 466 const Teuchos::ParameterList &pl = env->getParameters();
467 bool keep_partition_tree =
false;
468 const Teuchos::ParameterEntry * pe = pl.getEntryPtr(
"keep_partition_tree");
470 keep_partition_tree = pe->getValue(&keep_partition_tree);
471 if(!keep_partition_tree) {
472 throw std::logic_error(
473 "Requested tree when param keep_partition_tree is off.");
479 const Teuchos::ParameterList & zoltan_pl = pl.sublist(
"zoltan_parameters");
480 std::string lb_method;
481 pe = zoltan_pl.getEntryPtr(
"LB_METHOD");
483 lb_method = pe->getValue(&lb_method);
485 if(lb_method ==
"phg") {
486 phg_getPartitionTree(numParts, numTreeVerts, permPartNums,
487 splitRangeBeg, splitRangeEnd, treeVertParents);
489 else if(lb_method ==
"rcb") {
490 rcb_getPartitionTree(numParts, numTreeVerts, permPartNums,
491 splitRangeBeg, splitRangeEnd, treeVertParents);
494 throw std::logic_error(
"Did not recognize LB_METHOD: " + lb_method);
506 const RCP<
const Comm<int> > &problemComm__,
508 env(env__), problemComm(problemComm__), adapter(adapter__)
510 setMPIComm(problemComm__);
512 zz = rcp(
new Zoltan(mpicomm));
517 const RCP<
const Comm<int> > &problemComm__,
519 env(env__), problemComm(problemComm__), adapter(adapter__)
521 setMPIComm(problemComm__);
523 zz = rcp(
new Zoltan(mpicomm));
525 setCallbacksGeom(&(*adapter));
529 const RCP<
const Comm<int> > &problemComm__,
531 env(env__), problemComm(problemComm__), adapter(adapter__)
533 setMPIComm(problemComm__);
535 zz = rcp(
new Zoltan(mpicomm));
537 setCallbacksGraph(adapter);
538 if (adapter->coordinatesAvailable()) {
539 setCallbacksGeom(adapter->getCoordinateInput());
544 const RCP<
const Comm<int> > &problemComm__,
546 env(env__), problemComm(problemComm__), adapter(adapter__)
548 setMPIComm(problemComm__);
550 zz = rcp(
new Zoltan(mpicomm));
552 setCallbacksGraph(adapter);
553 setCallbacksHypergraph(adapter);
554 if (adapter->coordinatesAvailable()) {
555 setCallbacksGeom(adapter->getCoordinateInput());
560 const RCP<
const Comm<int> > &problemComm__,
562 env(env__), problemComm(problemComm__), adapter(adapter__)
564 setMPIComm(problemComm__);
566 zz = rcp(
new Zoltan(mpicomm));
568 setCallbacksGraph(adapter);
572 setCallbacksHypergraph(adapter);
573 setCallbacksGeom(&(*adapter));
581 template <
typename Adapter>
589 size_t numGlobalParts = solution->getTargetGlobalNumberOfParts();
591 sprintf(paramstr,
"%lu", numGlobalParts);
592 zz->Set_Param(
"NUM_GLOBAL_PARTS", paramstr);
594 int wdim = adapter->getNumWeightsPerID();
595 sprintf(paramstr,
"%d", wdim);
596 zz->Set_Param(
"OBJ_WEIGHT_DIM", paramstr);
598 const Teuchos::ParameterList &pl = env->getParameters();
601 const Teuchos::ParameterEntry *pe = pl.getEntryPtr(
"imbalance_tolerance");
604 tolerance = pe->getValue<
double>(&tolerance);
605 sprintf(str,
"%f", tolerance);
606 zz->Set_Param(
"IMBALANCE_TOL", str);
609 pe = pl.getEntryPtr(
"partitioning_approach");
611 std::string approach;
612 approach = pe->getValue<std::string>(&approach);
613 if (approach ==
"partition")
614 zz->Set_Param(
"LB_APPROACH",
"PARTITION");
616 zz->Set_Param(
"LB_APPROACH",
"REPARTITION");
619 pe = pl.getEntryPtr(
"partitioning_objective");
621 std::string strChoice = pe->getValue<std::string>(&strChoice);
622 if (strChoice == std::string(
"multicriteria_minimize_total_weight"))
623 zz->Set_Param(
"RCB_MULTICRITERIA_NORM",
"1");
624 else if (strChoice == std::string(
"multicriteria_balance_total_maximum"))
625 zz->Set_Param(
"RCB_MULTICRITERIA_NORM",
"2");
626 else if (strChoice == std::string(
"multicriteria_minimize_maximum_weight"))
627 zz->Set_Param(
"RCB_MULTICRITERIA_NORM",
"3");
634 bool keep_partition_tree =
false;
635 pe = pl.getEntryPtr(
"keep_partition_tree");
637 keep_partition_tree = pe->getValue(&keep_partition_tree);
638 if (keep_partition_tree) {
642 zz->Set_Param(
"KEEP_CUTS",
"1");
643 zz->Set_Param(
"PHG_KEEP_TREE",
"1");
647 pe = pl.getEntryPtr(
"rectilinear");
649 bool val = pe->getValue(&val);
651 zz->Set_Param(
"RCB_RECTILINEAR_BLOCKS",
"1");
656 const Teuchos::ParameterList &zpl = pl.sublist(
"zoltan_parameters");
657 for (ParameterList::ConstIterator iter = zpl.begin();
658 iter != zpl.end(); iter++) {
659 const std::string &zname = pl.name(iter);
661 std::string zval = pl.entry(iter).getValue(&zval);
662 zz->Set_Param(zname.c_str(), zval.c_str());
665 catch (std::exception &e) {
670 int pdim = (wdim > 1 ? wdim : 1);
671 for (
int d = 0; d < pdim; d++) {
672 if (!solution->criteriaHasUniformPartSizes(d)) {
673 float *partsizes =
new float[numGlobalParts];
674 int *partidx =
new int[numGlobalParts];
675 int *wgtidx =
new int[numGlobalParts];
676 for (
size_t i=0; i<numGlobalParts; i++) partidx[i] = i;
677 for (
size_t i=0; i<numGlobalParts; i++) wgtidx[i] = d;
678 for (
size_t i=0; i<numGlobalParts; i++)
679 partsizes[i] = solution->getCriteriaPartSize(0, i);
680 zz->LB_Set_Part_Sizes(1, numGlobalParts, partidx, wgtidx, partsizes);
693 ZOLTAN_ID_PTR dGids = NULL, dLids = NULL;
694 int *dProcs = NULL, *dParts = NULL;
696 ZOLTAN_ID_PTR oGids = NULL, oLids = NULL;
697 int *oProcs = NULL, *oParts = NULL;
699 zz->Set_Param(
"RETURN_LISTS",
"PARTS");
702 int ierr = zz->LB_Partition(changed, nGidEnt, nLidEnt,
703 nDummy, dGids, dLids, dProcs, dParts,
704 nObj, oGids, oLids, oProcs, oParts);
706 env->globalInputAssertion(__FILE__, __LINE__,
"Zoltan LB_Partition",
715 numObjects=model->getLocalNumObjects();
719 ArrayRCP<part_t> partList(
new part_t[numObjects], 0, numObjects,
true);
720 for (
int i = 0; i < nObj; i++) {
723 partList[tmp] = oParts[i];
735 Teuchos::RCP<const map_t> mapWithCopies;
736 Teuchos::RCP<const map_t> oneToOneMap;
739 typedef Tpetra::Vector<scalar_t, lno_t, gno_t> vector_t;
740 vector_t vecWithCopies(mapWithCopies);
741 vector_t oneToOneVec(oneToOneMap);
744 assert(nObj == lno_t(oneToOneMap->getNodeNumElements()));
745 for (lno_t i = 0; i < nObj; i++)
746 oneToOneVec.replaceLocalValue(i, oParts[i]);
749 Teuchos::RCP<const Tpetra::Import<lno_t, gno_t> > importer =
750 Tpetra::createImport<lno_t, gno_t>(oneToOneMap, mapWithCopies);
751 vecWithCopies.doImport(oneToOneVec, *importer, Tpetra::REPLACE);
754 lno_t nlocal = lno_t(mapWithCopies->getNodeNumElements());
755 for (lno_t i = 0; i < nlocal; i++)
756 partList[i] = vecWithCopies.getData()[i];
759 solution->setParts(partList);
762 zz->LB_Free_Part(&oGids, &oLids, &oProcs, &oParts);
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const VectorAdapter< user_t > > &adapter__)
IdentifierAdapter defines the interface for identifiers.
fast typical checks for valid arguments
MatrixAdapter defines the adapter interface for matrices.
Defines the Model interface.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const MatrixAdapter< user_t, userCoord_t > > &adapter__)
GraphAdapter defines the interface for graph-based user data.
MeshAdapter defines the interface for mesh input.
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
Defines the PartitioningSolution class.
SparseMatrixAdapter_t::part_t part_t
void partition(const RCP< PartitioningSolution< Adapter > > &solution)
Partitioning method.
callback functions for the Zoltan package (templated on Adapter)
A PartitioningSolution is a solution to a partitioning problem.
VectorAdapter defines the interface for vector input.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const MeshAdapter< user_t > > &adapter__)
Algorithm defines the base class for all algorithms.
static void ASSIGN(first_t &a, second_t b)
Traits class to handle conversions between gno_t/lno_t and TPL data types (e.g., ParMETIS's idx_t...
Gathering definitions used in software development.
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...
The base class for all model classes.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const GraphAdapter< user_t, userCoord_t > > &adapter__)
A gathering of useful namespace methods.
HyperGraphModel defines the interface required for hyper graph models.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const IdentifierAdapter< user_t > > &adapter__)