50 #ifndef _ZOLTAN2_PARTITIONINGPROBLEM_HPP_ 51 #define _ZOLTAN2_PARTITIONINGPROBLEM_HPP_ 62 #ifdef ZOLTAN2_TASKMAPPING_MOVE 74 #ifdef HAVE_ZOLTAN2_OVIS 103 template<
typename Adapter>
109 typedef typename Adapter::gno_t
gno_t;
110 typedef typename Adapter::lno_t
lno_t;
115 #ifdef HAVE_ZOLTAN2_MPI 116 typedef Teuchos::OpaqueWrapper<MPI_Comm> mpiWrapper_t;
120 Problem<Adapter>(A,p,comm), solution_(),
122 graphFlags_(), idFlags_(), coordFlags_(), algName_(),
123 numberOfWeights_(), partIds_(), partSizes_(),
124 numberOfCriteria_(), levelNumberParts_(), hierarchical_(false)
133 Problem<Adapter>(A,p), solution_(),
135 graphFlags_(), idFlags_(), coordFlags_(), algName_(),
137 partIds_(), partSizes_(), numberOfCriteria_(),
138 levelNumberParts_(), hierarchical_(false)
146 const RCP<
const Teuchos::Comm<int> > &comm):
147 Problem<Adapter>(A,p,comm), solution_(),
149 graphFlags_(), idFlags_(), coordFlags_(), algName_(),
151 partIds_(), partSizes_(), numberOfCriteria_(),
152 levelNumberParts_(), hierarchical_(false)
179 void solve(
bool updateInputData=
true );
186 return *(solution_.getRawPtr());
268 scalar_t *partSizes,
bool makeCopy=
true) ;
288 Array<std::string> algorithm_names(17);
289 algorithm_names[0] =
"rcb";
290 algorithm_names[1] =
"multijagged";
291 algorithm_names[2] =
"rib";
292 algorithm_names[3] =
"hsfc";
293 algorithm_names[4] =
"patoh";
294 algorithm_names[5] =
"phg";
295 algorithm_names[6] =
"metis";
296 algorithm_names[7] =
"parmetis";
297 algorithm_names[8] =
"pulp";
298 algorithm_names[9] =
"parma";
299 algorithm_names[10] =
"scotch";
300 algorithm_names[11] =
"ptscotch";
301 algorithm_names[12] =
"block";
302 algorithm_names[13] =
"cyclic";
303 algorithm_names[14] =
"random";
304 algorithm_names[15] =
"zoltan";
305 algorithm_names[16] =
"forTestingOnly";
306 RCP<Teuchos::StringValidator> algorithm_Validator = Teuchos::rcp(
307 new Teuchos::StringValidator( algorithm_names ));
308 pl.set(
"algorithm",
"random",
"partitioning algorithm",
309 algorithm_Validator);
312 pl.set(
"rectilinear",
false,
"If true, then when a cut is made, all of the " 313 "dots located on the cut are moved to the same side of the cut. The " 314 "resulting regions are then rectilinear. The resulting load balance may " 315 "not be as good as when the group of dots is split by the cut. ",
318 RCP<Teuchos::StringValidator> partitioning_objective_Validator =
319 Teuchos::rcp(
new Teuchos::StringValidator(
320 Teuchos::tuple<std::string>(
"balance_object_count",
321 "balance_object_weight",
"multicriteria_minimize_total_weight",
322 "multicriteria_minimize_maximum_weight",
323 "multicriteria_balance_total_maximum",
"minimize_cut_edge_count",
324 "minimize_cut_edge_weight",
"minimize_neighboring_parts",
325 "minimize_boundary_vertices" )));
326 pl.set(
"partitioning_objective",
"balance_object_weight",
327 "objective of partitioning", partitioning_objective_Validator);
329 pl.set(
"imbalance_tolerance", 1.1,
"imbalance tolerance, ratio of " 333 RCP<Teuchos::EnhancedNumberValidator<int>> num_global_parts_Validator =
334 Teuchos::rcp(
new Teuchos::EnhancedNumberValidator<int>(
335 1, Teuchos::EnhancedNumberTraits<int>::max()) );
336 pl.set(
"num_global_parts", 1,
"global number of parts to compute " 337 "(0 means use the number of processes)", num_global_parts_Validator);
340 RCP<Teuchos::EnhancedNumberValidator<int>> num_local_parts_Validator =
341 Teuchos::rcp(
new Teuchos::EnhancedNumberValidator<int>(
342 0, Teuchos::EnhancedNumberTraits<int>::max()) );
343 pl.set(
"num_local_parts", 0,
"number of parts to compute for this " 344 "process (0 means one)", num_local_parts_Validator);
346 RCP<Teuchos::StringValidator> partitioning_approach_Validator =
347 Teuchos::rcp(
new Teuchos::StringValidator(
348 Teuchos::tuple<std::string>(
"partition",
"repartition",
349 "maximize_overlap" )));
350 pl.set(
"partitioning_approach",
"partition",
"Partition from scratch, " 351 "partition incrementally from current partition, of partition from " 352 "scratch but maximize overlap with the current partition",
353 partitioning_approach_Validator);
355 RCP<Teuchos::StringValidator> objects_to_partition_Validator =
356 Teuchos::rcp(
new Teuchos::StringValidator(
357 Teuchos::tuple<std::string>(
"matrix_rows",
"matrix_columns",
358 "matrix_nonzeros",
"mesh_elements",
"mesh_nodes",
"graph_edges",
359 "graph_vertices",
"coordinates",
"identifiers" )));
360 pl.set(
"objects_to_partition",
"graph_vertices",
"Objects to be partitioned",
361 objects_to_partition_Validator);
363 RCP<Teuchos::StringValidator> model_Validator = Teuchos::rcp(
364 new Teuchos::StringValidator(
365 Teuchos::tuple<std::string>(
"hypergraph",
"graph",
366 "geometry",
"ids" )));
367 pl.set(
"model",
"graph",
"This is a low level parameter. Normally the " 368 "library will choose a computational model based on the algorithm or " 369 "objective specified by the user.", model_Validator);
372 pl.set(
"remap_parts",
false,
"remap part numbers to minimize migration " 375 pl.set(
"mapping_type", -1,
"Mapping of solution to the processors. -1 No" 378 RCP<Teuchos::EnhancedNumberValidator<int>> ghost_layers_Validator =
379 Teuchos::rcp(
new Teuchos::EnhancedNumberValidator<int>(1, 10, 1, 0) );
380 pl.set(
"ghost_layers", 2,
"number of layers for ghosting used in " 381 "hypergraph ghost method", ghost_layers_Validator);
385 void initializeProblem();
387 void createPartitioningProblem(
bool newData);
389 RCP<PartitioningSolution<Adapter> > solution_;
390 #ifdef ZOLTAN2_TASKMAPPING_MOVE 391 RCP<MachineRepresentation<scalar_t,part_t> > machine_;
402 std::string algName_;
404 int numberOfWeights_;
415 ArrayRCP<ArrayRCP<part_t> > partIds_;
416 ArrayRCP<ArrayRCP<scalar_t> > partSizes_;
417 int numberOfCriteria_;
421 ArrayRCP<int> levelNumberParts_;
434 template <
typename Adapter>
441 if (getenv(
"DEBUGME")){
443 std::cout << getpid() << std::endl;
446 std::cout << _getpid() << std::endl;
451 #ifdef HAVE_ZOLTAN2_OVIS 452 ovis_enabled(this->
comm_->getRank());
457 #ifdef ZOLTAN2_TASKMAPPING_MOVE 458 machine_ = RCP<MachineRepresentation<scalar_t,part_t> >(
465 numberOfWeights_ = this->
inputAdapter_->getNumWeightsPerID();
467 numberOfCriteria_ = (numberOfWeights_ > 1) ? numberOfWeights_ : 1;
474 ArrayRCP<part_t> *noIds =
new ArrayRCP<part_t> [numberOfCriteria_];
475 ArrayRCP<scalar_t> *noSizes =
new ArrayRCP<scalar_t> [numberOfCriteria_];
477 partIds_ = arcp(noIds, 0, numberOfCriteria_,
true);
478 partSizes_ = arcp(noSizes, 0, numberOfCriteria_,
true);
481 std::ostringstream msg;
482 msg << this->
comm_->getSize() <<
" procs," 483 << numberOfWeights_ <<
" user-defined weights\n";
487 this->
env_->memory(
"After initializeProblem");
490 template <
typename Adapter>
492 int criteria,
int len,
part_t *partIds,
scalar_t *partSizes,
bool makeCopy)
494 this->
env_->localInputAssertion(__FILE__, __LINE__,
"invalid length",
497 this->
env_->localInputAssertion(__FILE__, __LINE__,
"invalid criteria",
501 partIds_[criteria] = ArrayRCP<part_t>();
502 partSizes_[criteria] = ArrayRCP<scalar_t>();
506 this->
env_->localInputAssertion(__FILE__, __LINE__,
"invalid arrays",
513 part_t *z2_partIds = NULL;
515 bool own_memory =
false;
518 z2_partIds =
new part_t [len];
520 this->
env_->localMemoryAssertion(__FILE__, __LINE__, len, z2_partSizes);
521 memcpy(z2_partIds, partIds, len *
sizeof(
part_t));
522 memcpy(z2_partSizes, partSizes, len *
sizeof(
scalar_t));
526 z2_partIds = partIds;
527 z2_partSizes = partSizes;
530 partIds_[criteria] = arcp(z2_partIds, 0, len, own_memory);
531 partSizes_[criteria] = arcp(z2_partSizes, 0, len, own_memory);
534 template <
typename Adapter>
543 createPartitioningProblem(updateInputData);
555 if (algName_ == std::string(
"multijagged")) {
560 else if (algName_ == std::string(
"zoltan")) {
565 else if (algName_ == std::string(
"parma")) {
570 else if (algName_ == std::string(
"scotch")) {
575 else if (algName_ == std::string(
"parmetis")) {
580 else if (algName_ == std::string(
"pulp")) {
585 else if (algName_ == std::string(
"block")) {
589 else if (algName_ == std::string(
"phg") ||
590 algName_ == std::string(
"patoh")) {
592 Teuchos::ParameterList &pl = this->
env_->getParametersNonConst();
593 Teuchos::ParameterList &zparams = pl.sublist(
"zoltan_parameters",
false);
594 if (numberOfWeights_ > 0) {
596 sprintf(strval,
"%d", numberOfWeights_);
597 zparams.set(
"OBJ_WEIGHT_DIM", strval);
599 zparams.set(
"LB_METHOD", algName_.c_str());
600 zparams.set(
"LB_APPROACH",
"PARTITION");
601 algName_ = std::string(
"zoltan");
607 else if (algName_ == std::string(
"forTestingOnly")) {
617 throw std::logic_error(
"partitioning algorithm not supported");
629 partIds_.view(0, numberOfCriteria_),
630 partSizes_.view(0, numberOfCriteria_), this->
algorithm_);
634 solution_ = rcp(soln);
637 this->
env_->memory(
"After creating Solution");
647 const Teuchos::ParameterEntry *pe = this->
envConst_->getParameters().getEntryPtr(
"mapping_type");
648 int mapping_type = -1;
650 mapping_type = pe->getValue(&mapping_type);
654 #if ZOLTAN2_TASKMAPPING_MOVE 655 if (mapping_type == 0){
661 this->
comm_.getRawPtr(),
662 machine_.getRawPtr(),
664 solution_.getRawPtr(),
678 const part_t *oldParts = solution_->getPartListView();
679 size_t nLocal = ia->getNumLocalIds();
680 for (
size_t i = 0; i < nLocal; i++) {
692 else if (mapping_type == 1){
699 template <
typename Adapter>
703 "PartitioningProblem::createPartitioningProblem");
705 using Teuchos::ParameterList;
716 prevModelAvail[i] = modelAvail_[i];
721 modelFlag_t previousIdentifierModelFlags = idFlags_;
722 modelFlag_t previousCoordinateModelFlags = coordFlags_;
727 modelAvail_[i] =
false;
757 std::string defString(
"default");
761 std::string model(defString);
762 const Teuchos::ParameterEntry *pe = pl.getEntryPtr(
"model");
764 model = pe->getValue<std::string>(&model);
768 std::string algorithm(defString);
769 pe = pl.getEntryPtr(
"algorithm");
771 algorithm = pe->getValue<std::string>(&algorithm);
775 bool needConsecutiveGlobalIds =
false;
776 bool removeSelfEdges=
false;
782 if (algorithm != defString)
786 if (algorithm == std::string(
"block") ||
787 algorithm == std::string(
"random") ||
788 algorithm == std::string(
"cyclic") ){
793 algName_ = algorithm;
795 else if (algorithm == std::string(
"zoltan") ||
796 algorithm == std::string(
"parma") ||
797 algorithm == std::string(
"forTestingOnly"))
799 algName_ = algorithm;
801 else if (algorithm == std::string(
"rcb") ||
802 algorithm == std::string(
"rib") ||
803 algorithm == std::string(
"hsfc"))
806 Teuchos::ParameterList &zparams = pl.sublist(
"zoltan_parameters",
false);
807 zparams.set(
"LB_METHOD", algorithm);
808 if (numberOfWeights_ > 0) {
810 sprintf(strval,
"%d", numberOfWeights_);
811 zparams.set(
"OBJ_WEIGHT_DIM", strval);
813 algName_ = std::string(
"zoltan");
815 else if (algorithm == std::string(
"multijagged"))
820 algName_ = algorithm;
822 else if (algorithm == std::string(
"metis") ||
823 algorithm == std::string(
"parmetis"))
828 algName_ = algorithm;
829 removeSelfEdges =
true;
830 needConsecutiveGlobalIds =
true;
832 else if (algorithm == std::string(
"scotch") ||
833 algorithm == std::string(
"ptscotch"))
835 algName_ = algorithm;
837 else if (algorithm == std::string(
"pulp"))
839 algName_ = algorithm;
841 else if (algorithm == std::string(
"patoh") ||
842 algorithm == std::string(
"phg"))
852 algName_ = algorithm;
854 #ifdef INCLUDE_ZOLTAN2_EXPERIMENTAL_WOLF 855 else if (algorithm == std::string(
"nd"))
859 algName_ = algorithm;
865 throw std::logic_error(
"parameter list algorithm is invalid");
868 else if (model != defString)
871 if (model == std::string(
"hypergraph"))
876 if (this->
comm_->getSize() > 1)
877 algName_ = std::string(
"phg");
879 algName_ = std::string(
"patoh");
881 else if (model == std::string(
"graph"))
886 #ifdef HAVE_ZOLTAN2_SCOTCH 888 if (this->
comm_->getSize() > 1)
889 algName_ = std::string(
"ptscotch");
891 algName_ = std::string(
"scotch");
893 #ifdef HAVE_ZOLTAN2_PARMETIS 894 if (this->
comm_->getSize() > 1)
895 algName_ = std::string(
"parmetis");
897 algName_ = std::string(
"metis");
898 removeSelfEdges =
true;
899 needConsecutiveGlobalIds =
true;
901 #ifdef HAVE_ZOLTAN2_PULP 906 algName_ = std::string(
"pulp");
908 if (this->
comm_->getSize() > 1)
909 algName_ = std::string(
"phg");
911 algName_ = std::string(
"patoh");
916 else if (model == std::string(
"geometry"))
921 algName_ = std::string(
"multijagged");
923 else if (model == std::string(
"ids"))
928 algName_ = std::string(
"block");
948 if (this->
comm_->getSize() > 1)
949 algName_ = std::string(
"phg");
951 algName_ = std::string(
"patoh");
959 if (this->
comm_->getSize() > 1)
960 algName_ = std::string(
"phg");
962 algName_ = std::string(
"patoh");
969 algName_ = std::string(
"multijagged");
976 algName_ = std::string(
"block");
980 throw std::logic_error(
"input type is invalid");
986 Array<int> valueList;
987 pe = pl.getEntryPtr(
"topology");
990 valueList = pe->getValue<Array<int> >(&valueList);
992 if (!Zoltan2::noValuesAreInRangeList<int>(valueList)){
993 int *n =
new int [valueList.size() + 1];
994 levelNumberParts_ = arcp(n, 0, valueList.size() + 1,
true);
995 int procsPerNode = 1;
996 for (
int i=0; i < valueList.size(); i++){
997 levelNumberParts_[i+1] = valueList[i];
998 procsPerNode *= valueList[i];
1001 levelNumberParts_[0] = env.
numProcs_ / procsPerNode;
1004 levelNumberParts_[0]++;
1008 levelNumberParts_.clear();
1011 hierarchical_ = levelNumberParts_.size() > 0;
1015 std::string objectOfInterest(defString);
1016 pe = pl.getEntryPtr(
"objects_to_partition");
1018 objectOfInterest = pe->getValue<std::string>(&objectOfInterest);
1030 std::string symParameter(defString);
1031 pe = pl.getEntryPtr(
"symmetrize_graph");
1033 symParameter = pe->getValue<std::string>(&symParameter);
1034 if (symParameter == std::string(
"transpose"))
1036 else if (symParameter == std::string(
"bipartite"))
1040 bool sgParameter =
false;
1041 pe = pl.getEntryPtr(
"subset_graph");
1043 sgParameter = pe->getValue(&sgParameter);
1045 if (sgParameter == 1)
1050 if (removeSelfEdges)
1053 if (needConsecutiveGlobalIds)
1059 if (objectOfInterest == defString ||
1060 objectOfInterest == std::string(
"matrix_rows") )
1062 else if (objectOfInterest == std::string(
"matrix_columns"))
1064 else if (objectOfInterest == std::string(
"matrix_nonzeros"))
1069 if (objectOfInterest == defString ||
1070 objectOfInterest == std::string(
"mesh_nodes") )
1072 else if (objectOfInterest == std::string(
"mesh_elements"))
1099 (graphFlags_ != previousGraphModelFlags) ||
1100 (coordFlags_ != previousCoordinateModelFlags) ||
1101 (idFlags_ != previousIdentifierModelFlags))
1113 if(modelAvail_[GraphModelType]==
true)
1123 if(modelAvail_[HypergraphModelType]==
true)
1129 if(modelAvail_[CoordinateModelType]==
true)
1140 if(modelAvail_[IdentifierModelType]==
true)
1151 this->
env_->memory(
"After creating Model");
use columns as graph vertices
Serial greedy first-fit graph coloring (local graph only)
RCP< GraphModel< base_adapter_t > > graphModel_
~PartitioningProblem()
Destructor.
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Time an algorithm (or other entity) as a whole.
RCP< const base_adapter_t > baseInputAdapter_
Adapter::scalar_t scalar_t
void setPartSizes(int len, part_t *partIds, scalar_t *partSizes, bool makeCopy=true)
Set or reset relative sizes for the parts that Zoltan2 will create.
use mesh nodes as vertices
fast typical checks for valid arguments
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
static void getValidParameters(ParameterList &pl)
Set up validators specific to this algorithm.
static RCP< Teuchos::BoolParameterEntryValidator > getBoolValidator()
Exists to make setting up validators less cluttered.
algorithm requires consecutive ids
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
PartitioningProblem(Adapter *A, ParameterList *p)
Constructor where communicator is the Teuchos default.
model must symmetrize input
model must symmetrize input
RCP< Algorithm< Adapter > > algorithm_
PartitioningProblem(Adapter *A, ParameterList *p, const RCP< const Teuchos::Comm< int > > &comm)
Constructor where Teuchos communicator is specified.
RCP< const Comm< int > > comm_
Defines the PartitioningSolution class.
static RCP< Teuchos::AnyNumberParameterEntryValidator > getAnyDoubleValidator()
Exists to make setting up validators less cluttered.
part_t getAssignedProcForTask(part_t taskId)
getAssignedProcForTask function, returns the assigned processor id for the given task ...
sub-steps, each method's entry and exit
static void getValidParameters(ParameterList &pl)
Set up validators specific to this Problem.
RCP< IdentifierModel< base_adapter_t > > identifierModel_
void setPartSizesForCriteria(int criteria, int len, part_t *partIds, scalar_t *partSizes, bool makeCopy=true)
Set or reset the relative sizes (per weight) for the parts that Zoltan2 will create.
int numProcs_
number of processes (relative to comm_)
use matrix rows as graph vertices
This class provides geometric coordinates with optional weights to the Zoltan2 algorithm.
Teuchos::ParameterList & getParametersNonConst()
Returns a reference to a non-const copy of the parameters.
algorithm requires no self edges
Defines the IdentifierModel interface.
A PartitioningSolution is a solution to a partitioning problem.
use nonzeros as graph vertices
Defines the EvaluatePartition class.
void localBugAssertion(const char *file, int lineNum, const char *msg, bool ok, AssertionLevel level) const
Test for valid library behavior on local process only.
BaseAdapterType
An enum to identify general types of adapters.
identifier data, just a list of IDs
Problem base class from which other classes (PartitioningProblem, ColoringProblem, OrderingProblem, MatchingProblem, etc.) derive.
RCP< const Adapter > inputAdapter_
Defines the Problem base class.
RCP< CoordinateModel< base_adapter_t > > coordinateModel_
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
static RCP< Teuchos::AnyNumberParameterEntryValidator > getAnyIntValidator()
Exists to make setting up validators less cluttered.
static void getValidParameters(ParameterList &pl)
Set up validators specific to this algorithm.
MachineRepresentation Class Base class for representing machine coordinates, networks, etc.
Adapter::base_adapter_t base_adapter_t
Define IntegerRangeList validator.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
static void getValidParameters(ParameterList &pl)
Set up validators specific to this algorithm.
PartitioningProblem sets up partitioning problems for the user.
use mesh elements as vertices
GraphModel defines the interface required for graph models.
The base class for all model classes.
IdentifierModel defines the interface for all identifier models.
Defines the GraphModel interface.
static void getValidParameters(ParameterList &pl)
Set up validators specific to this algorithm.
RCP< const Model< base_adapter_t > > baseModel_
void solve(bool updateInputData=true)
Direct the problem to create a solution.
RCP< const Environment > envConst_
Multi Jagged coordinate partitioning algorithm.
static void getValidParameters(ParameterList &pl)
Set up validators specific to this algorithm.