Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_PartitioningProblem.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
50 #ifndef _ZOLTAN2_PARTITIONINGPROBLEM_HPP_
51 #define _ZOLTAN2_PARTITIONINGPROBLEM_HPP_
52 
53 #include <Zoltan2_Problem.hpp>
57 #include <Zoltan2_GraphModel.hpp>
62 #ifdef ZOLTAN2_TASKMAPPING_MOVE
63 #include <Zoltan2_TaskMapping.hpp>
64 #endif
65 
66 #ifndef _WIN32
67 #include <unistd.h>
68 #else
69 #include <process.h>
70 #define NOMINMAX
71 #include <windows.h>
72 #endif
73 
74 #ifdef HAVE_ZOLTAN2_OVIS
75 #include <ovis.h>
76 #endif
77 
78 namespace Zoltan2{
79 
103 template<typename Adapter>
104 class PartitioningProblem : public Problem<Adapter>
105 {
106 public:
107 
108  typedef typename Adapter::scalar_t scalar_t;
109  typedef typename Adapter::gno_t gno_t;
110  typedef typename Adapter::lno_t lno_t;
111  typedef typename Adapter::part_t part_t;
112  typedef typename Adapter::user_t user_t;
114 
116  PartitioningProblem(Adapter *A, ParameterList *p,
117  const RCP<const Teuchos::Comm<int> > &comm):
118  Problem<Adapter>(A,p,comm),
119  solution_(),
120  inputType_(InvalidAdapterType),
121  graphFlags_(), idFlags_(), coordFlags_(),
122  algName_(), numberOfWeights_(), partIds_(), partSizes_(),
123  numberOfCriteria_(), levelNumberParts_(), hierarchical_(false)
124  {
125  for(int i=0;i<MAX_NUM_MODEL_TYPES;i++) modelAvail_[i]=false;
126  initializeProblem();
127  }
128 
129 #ifdef HAVE_ZOLTAN2_MPI
130 
132  PartitioningProblem(Adapter *A, ParameterList *p, MPI_Comm mpicomm):
133  PartitioningProblem(A, p,
134  rcp<const Comm<int> >(new Teuchos::MpiComm<int>(
135  Teuchos::opaqueWrapper(mpicomm))))
136  {}
137 #endif
138 
140  PartitioningProblem(Adapter *A, ParameterList *p):
141  PartitioningProblem(A, p, Tpetra::getDefaultComm())
142  {}
143 
147 
149  //
150  // \param updateInputData If true this indicates that either
151  // this is the first attempt at solution, or that we
152  // are computing a new solution and the input data has
153  // changed since the previous solution was computed.
154  // By input data we mean coordinates, topology, or weights.
155  // If false, this indicates that we are computing a
156  // new solution using the same input data was used for
157  // the previous solution, even though the parameters
158  // may have been changed.
159  //
160  // For the sake of performance, we ask the caller to set \c updateInputData
161  // to false if he/she is computing a new solution using the same input data,
162  // but different problem parameters, than that which was used to compute
163  // the most recent solution.
164 
165  void solve(bool updateInputData=true);
166 
168  //
169  // \return a reference to the solution to the most recent solve().
170 
172  return *(solution_.getRawPtr());
173  };
174 
213  void setPartSizes(int len, part_t *partIds, scalar_t *partSizes,
214  bool makeCopy=true)
215  {
216  setPartSizesForCriteria(0, len, partIds, partSizes, makeCopy);
217  }
218 
253  void setPartSizesForCriteria(int criteria, int len, part_t *partIds,
254  scalar_t *partSizes, bool makeCopy=true) ;
255 /*
256  void setMachine(MachineRepresentation<typename Adapter::base_adapter_t::scalar_t> *machine);
257 */
258 
261  static void getValidParameters(ParameterList & pl)
262  {
268 
269  // This set up does not use tuple because we didn't have constructors
270  // that took that many elements - Tuple will need to be modified and I
271  // didn't want to have low level changes with this particular refactor
272  // TO DO: Add more Tuple constructors and then redo this code to be
273  // Teuchos::tuple<std::string> algorithm_names( "rcb", "multijagged" ... );
274  Array<std::string> algorithm_names(17);
275  algorithm_names[0] = "rcb";
276  algorithm_names[1] = "multijagged";
277  algorithm_names[2] = "rib";
278  algorithm_names[3] = "hsfc";
279  algorithm_names[4] = "patoh";
280  algorithm_names[5] = "phg";
281  algorithm_names[6] = "metis";
282  algorithm_names[7] = "parmetis";
283  algorithm_names[8] = "pulp";
284  algorithm_names[9] = "parma";
285  algorithm_names[10] = "scotch";
286  algorithm_names[11] = "ptscotch";
287  algorithm_names[12] = "block";
288  algorithm_names[13] = "cyclic";
289  algorithm_names[14] = "random";
290  algorithm_names[15] = "zoltan";
291  algorithm_names[16] = "forTestingOnly";
292  RCP<Teuchos::StringValidator> algorithm_Validator = Teuchos::rcp(
293  new Teuchos::StringValidator( algorithm_names ));
294  pl.set("algorithm", "random", "partitioning algorithm",
295  algorithm_Validator);
296 
297  // bool parameter
298  pl.set("keep_partition_tree", false, "If true, will keep partition tree",
300 
301  // bool parameter
302  pl.set("rectilinear", false, "If true, then when a cut is made, all of the "
303  "dots located on the cut are moved to the same side of the cut. The "
304  "resulting regions are then rectilinear. The resulting load balance may "
305  "not be as good as when the group of dots is split by the cut. ",
307 
308  RCP<Teuchos::StringValidator> partitioning_objective_Validator =
309  Teuchos::rcp( new Teuchos::StringValidator(
310  Teuchos::tuple<std::string>( "balance_object_count",
311  "balance_object_weight", "multicriteria_minimize_total_weight",
312  "multicriteria_minimize_maximum_weight",
313  "multicriteria_balance_total_maximum", "minimize_cut_edge_count",
314  "minimize_cut_edge_weight", "minimize_neighboring_parts",
315  "minimize_boundary_vertices" )));
316  pl.set("partitioning_objective", "balance_object_weight",
317  "objective of partitioning", partitioning_objective_Validator);
318 
319  pl.set("imbalance_tolerance", 1.1, "imbalance tolerance, ratio of "
320  "maximum load over average load", Environment::getAnyDoubleValidator());
321 
322  // num_global_parts >= 1
323  RCP<Teuchos::EnhancedNumberValidator<int>> num_global_parts_Validator =
324  Teuchos::rcp( new Teuchos::EnhancedNumberValidator<int>(
325  1, Teuchos::EnhancedNumberTraits<int>::max()) ); // no maximum
326  pl.set("num_global_parts", 1, "global number of parts to compute "
327  "(0 means use the number of processes)", num_global_parts_Validator);
328 
329  // num_local_parts >= 0
330  RCP<Teuchos::EnhancedNumberValidator<int>> num_local_parts_Validator =
331  Teuchos::rcp( new Teuchos::EnhancedNumberValidator<int>(
332  0, Teuchos::EnhancedNumberTraits<int>::max()) ); // no maximum
333  pl.set("num_local_parts", 0, "number of parts to compute for this "
334  "process (num_global_parts == sum of all num_local_parts)",
335  num_local_parts_Validator);
336 
337  RCP<Teuchos::StringValidator> partitioning_approach_Validator =
338  Teuchos::rcp( new Teuchos::StringValidator(
339  Teuchos::tuple<std::string>( "partition", "repartition",
340  "maximize_overlap" )));
341  pl.set("partitioning_approach", "partition", "Partition from scratch, "
342  "partition incrementally from current partition, of partition from "
343  "scratch but maximize overlap with the current partition",
344  partitioning_approach_Validator);
345 
346  RCP<Teuchos::StringValidator> objects_to_partition_Validator =
347  Teuchos::rcp( new Teuchos::StringValidator(
348  Teuchos::tuple<std::string>( "matrix_rows", "matrix_columns",
349  "matrix_nonzeros", "mesh_elements", "mesh_nodes", "graph_edges",
350  "graph_vertices", "coordinates", "identifiers" )));
351  pl.set("objects_to_partition", "graph_vertices", "Objects to be partitioned",
352  objects_to_partition_Validator);
353 
354  RCP<Teuchos::StringValidator> model_Validator = Teuchos::rcp(
355  new Teuchos::StringValidator(
356  Teuchos::tuple<std::string>( "hypergraph", "graph",
357  "geometry", "ids" )));
358  pl.set("model", "graph", "This is a low level parameter. Normally the "
359  "library will choose a computational model based on the algorithm or "
360  "objective specified by the user.", model_Validator);
361 
362  // bool parameter
363  pl.set("remap_parts", false, "remap part numbers to minimize migration "
364  "between old and new partitions", Environment::getBoolValidator() );
365 
366  pl.set("mapping_type", -1, "Mapping of solution to the processors. -1 No"
367  " Mapping, 0 coordinate mapping.", Environment::getAnyIntValidator());
368 
369  RCP<Teuchos::EnhancedNumberValidator<int>> ghost_layers_Validator =
370  Teuchos::rcp( new Teuchos::EnhancedNumberValidator<int>(1, 10, 1, 0) );
371  pl.set("ghost_layers", 2, "number of layers for ghosting used in "
372  "hypergraph ghost method", ghost_layers_Validator);
373  }
374 
375 private:
376  void initializeProblem();
377 
378  void createPartitioningProblem(bool newData);
379 
380  RCP<PartitioningSolution<Adapter> > solution_;
381 #ifdef ZOLTAN2_TASKMAPPING_MOVE
382  RCP<MachineRepresentation<scalar_t,part_t> > machine_;
383 #endif
384 
385  BaseAdapterType inputType_;
386 
387  //ModelType modelType_;
388  bool modelAvail_[MAX_NUM_MODEL_TYPES];
389 
390  modelFlag_t graphFlags_;
391  modelFlag_t idFlags_;
392  modelFlag_t coordFlags_;
393  std::string algName_;
394 
395  int numberOfWeights_;
396 
397  // Suppose Array<part_t> partIds = partIds_[w]. If partIds.size() > 0
398  // then the user supplied part sizes for weight index "w", and the sizes
399  // corresponding to the Ids in partIds are partSizes[w].
400  //
401  // If numberOfWeights_ >= 0, then there is an Id and Sizes array for
402  // for each weight. Otherwise the user did not supply object weights,
403  // but they can still specify part sizes.
404  // So numberOfCriteria_ is numberOfWeights_ or one, whichever is greater.
405 
406  ArrayRCP<ArrayRCP<part_t> > partIds_;
407  ArrayRCP<ArrayRCP<scalar_t> > partSizes_;
408  int numberOfCriteria_;
409 
410  // Number of parts to be computed at each level in hierarchical partitioning.
411 
412  ArrayRCP<int> levelNumberParts_;
413  bool hierarchical_;
414 };
416 
417 /*
418 template <typename Adapter>
419 void PartitioningProblem<Adapter>::setMachine(MachineRepresentation<typename Adapter::base_adapter_t::scalar_t> *machine){
420  this->machine_ = RCP<MachineRepresentation<typename Adapter::base_adapter_t::scalar_t> > (machine, false);
421 }
422 */
423 
424 
425 template <typename Adapter>
426  void PartitioningProblem<Adapter>::initializeProblem()
427 {
428  HELLO;
429 
430  this->env_->debug(DETAILED_STATUS, "PartitioningProblem::initializeProblem");
431 
432  if (getenv("DEBUGME")){
433 #ifndef _WIN32
434  std::cout << getpid() << std::endl;
435  sleep(15);
436 #else
437  std::cout << _getpid() << std::endl;
438  Sleep(15000);
439 #endif
440  }
441 
442 #ifdef HAVE_ZOLTAN2_OVIS
443  ovis_enabled(this->comm_->getRank());
444 #endif
445 
446  // Create a copy of the user's communicator.
447 
448 #ifdef ZOLTAN2_TASKMAPPING_MOVE
449  machine_ = RCP<MachineRepresentation<scalar_t,part_t> >(
450  new MachineRepresentation<scalar_t,part_t>(*(this->comm_)));
451 #endif
452 
453  // Number of criteria is number of user supplied weights if non-zero.
454  // Otherwise it is 1 and uniform weight is implied.
455 
456  numberOfWeights_ = this->inputAdapter_->getNumWeightsPerID();
457 
458  numberOfCriteria_ = (numberOfWeights_ > 1) ? numberOfWeights_ : 1;
459 
460  inputType_ = this->inputAdapter_->adapterType();
461 
462  // The Caller can specify part sizes in setPartSizes(). If he/she
463  // does not, the part size arrays are empty.
464 
465  ArrayRCP<part_t> *noIds = new ArrayRCP<part_t> [numberOfCriteria_];
466  ArrayRCP<scalar_t> *noSizes = new ArrayRCP<scalar_t> [numberOfCriteria_];
467 
468  partIds_ = arcp(noIds, 0, numberOfCriteria_, true);
469  partSizes_ = arcp(noSizes, 0, numberOfCriteria_, true);
470 
471  if (this->env_->getDebugLevel() >= DETAILED_STATUS){
472  std::ostringstream msg;
473  msg << this->comm_->getSize() << " procs,"
474  << numberOfWeights_ << " user-defined weights\n";
475  this->env_->debug(DETAILED_STATUS, msg.str());
476  }
477 
478  this->env_->memory("After initializeProblem");
479 }
480 
481 template <typename Adapter>
483  int criteria, int len, part_t *partIds, scalar_t *partSizes, bool makeCopy)
484 {
485  this->env_->localInputAssertion(__FILE__, __LINE__, "invalid length",
486  len>= 0, BASIC_ASSERTION);
487 
488  this->env_->localInputAssertion(__FILE__, __LINE__, "invalid criteria",
489  criteria >= 0 && criteria < numberOfCriteria_, BASIC_ASSERTION);
490 
491  if (len == 0){
492  partIds_[criteria] = ArrayRCP<part_t>();
493  partSizes_[criteria] = ArrayRCP<scalar_t>();
494  return;
495  }
496 
497  this->env_->localInputAssertion(__FILE__, __LINE__, "invalid arrays",
498  partIds && partSizes, BASIC_ASSERTION);
499 
500  // The global validity of the partIds and partSizes arrays is performed
501  // by the PartitioningSolution, which computes global part distribution and
502  // part sizes.
503 
504  part_t *z2_partIds = NULL;
505  scalar_t *z2_partSizes = NULL;
506  bool own_memory = false;
507 
508  if (makeCopy){
509  z2_partIds = new part_t [len];
510  z2_partSizes = new scalar_t [len];
511  this->env_->localMemoryAssertion(__FILE__, __LINE__, len, z2_partSizes);
512  memcpy(z2_partIds, partIds, len * sizeof(part_t));
513  memcpy(z2_partSizes, partSizes, len * sizeof(scalar_t));
514  own_memory=true;
515  }
516  else{
517  z2_partIds = partIds;
518  z2_partSizes = partSizes;
519  }
520 
521  partIds_[criteria] = arcp(z2_partIds, 0, len, own_memory);
522  partSizes_[criteria] = arcp(z2_partSizes, 0, len, own_memory);
523 }
524 
525 template <typename Adapter>
526 void PartitioningProblem<Adapter>::solve(bool updateInputData)
527 {
528  this->env_->debug(DETAILED_STATUS, "Entering solve");
529 
530  // Create the computational model.
531 
532  this->env_->timerStart(MACRO_TIMERS, "create problem");
533 
534  createPartitioningProblem(updateInputData);
535 
536  this->env_->timerStop(MACRO_TIMERS, "create problem");
537 
538  // TODO: If hierarchical_
539 
540  // Create the solution. The algorithm will query the Solution
541  // for part and weight information. The algorithm will
542  // update the solution with part assignments and quality metrics.
543 
544  // Create the algorithm
545  try {
546  if (algName_ == std::string("multijagged")) {
547  this->algorithm_ = rcp(new Zoltan2_AlgMJ<Adapter>(this->envConst_,
548  this->comm_,
549  this->coordinateModel_));
550  }
551  else if (algName_ == std::string("zoltan")) {
552  this->algorithm_ = rcp(new AlgZoltan<Adapter>(this->envConst_,
553  this->comm_,
554  this->baseInputAdapter_));
555  }
556  else if (algName_ == std::string("parma")) {
557  this->algorithm_ = rcp(new AlgParMA<Adapter>(this->envConst_,
558  this->comm_,
559  this->baseInputAdapter_));
560  }
561  else if (algName_ == std::string("scotch")) {
562  this->algorithm_ = rcp(new AlgPTScotch<Adapter>(this->envConst_,
563  this->comm_,
564  this->baseInputAdapter_));
565  }
566  else if (algName_ == std::string("parmetis")) {
567  this->algorithm_ = rcp(new AlgParMETIS<Adapter>(this->envConst_,
568  this->comm_,
569  this->graphModel_));
570  }
571  else if (algName_ == std::string("pulp")) {
572  this->algorithm_ = rcp(new AlgPuLP<Adapter>(this->envConst_,
573  this->comm_,
574  this->baseInputAdapter_));
575  }
576  else if (algName_ == std::string("block")) {
577  this->algorithm_ = rcp(new AlgBlock<Adapter>(this->envConst_,
578  this->comm_, this->identifierModel_));
579  }
580  else if (algName_ == std::string("phg") ||
581  algName_ == std::string("patoh")) {
582  // phg and patoh provided through Zoltan
583  Teuchos::ParameterList &pl = this->env_->getParametersNonConst();
584  Teuchos::ParameterList &zparams = pl.sublist("zoltan_parameters",false);
585  if (numberOfWeights_ > 0) {
586  char strval[20];
587  sprintf(strval, "%d", numberOfWeights_);
588  zparams.set("OBJ_WEIGHT_DIM", strval);
589  }
590  zparams.set("LB_METHOD", algName_.c_str());
591  zparams.set("LB_APPROACH", "PARTITION");
592  algName_ = std::string("zoltan");
593 
594  this->algorithm_ = rcp(new AlgZoltan<Adapter>(this->envConst_,
595  this->comm_,
596  this->baseInputAdapter_));
597  }
598  else if (algName_ == std::string("forTestingOnly")) {
599  this->algorithm_ = rcp(new AlgForTestingOnly<Adapter>(this->envConst_,
600  this->comm_,
601  this->baseInputAdapter_));
602  }
603  // else if (algName_ == std::string("rcb")) {
604  // this->algorithm_ = rcp(new AlgRCB<Adapter>(this->envConst_,this->comm_,
605  // this->coordinateModel_));
606  // }
607  else {
608  throw std::logic_error("partitioning algorithm not supported");
609  }
610  }
612 
613  // Create the solution
614  this->env_->timerStart(MACRO_TIMERS, "create solution");
615  PartitioningSolution<Adapter> *soln = NULL;
616 
617  try{
619  this->envConst_, this->comm_, numberOfWeights_,
620  partIds_.view(0, numberOfCriteria_),
621  partSizes_.view(0, numberOfCriteria_), this->algorithm_);
622  }
624 
625  solution_ = rcp(soln);
626 
627  this->env_->timerStop(MACRO_TIMERS, "create solution");
628  this->env_->memory("After creating Solution");
629 
630  // Call the algorithm
631 
632  try {
633  this->algorithm_->partition(solution_);
634  }
636 
637  //if mapping is requested
638  const Teuchos::ParameterEntry *pe = this->envConst_->getParameters().getEntryPtr("mapping_type");
639  int mapping_type = -1;
640  if (pe){
641  mapping_type = pe->getValue(&mapping_type);
642  }
643  //if mapping is 0 -- coordinate mapping
644 
645 #if ZOLTAN2_TASKMAPPING_MOVE
646  if (mapping_type == 0){
647 
648  //part_t *task_communication_xadj = NULL, *task_communication_adj = NULL;
649 
652  this->comm_.getRawPtr(),
653  machine_.getRawPtr(),
654  this->coordinateModel_.getRawPtr(),
655  solution_.getRawPtr(),
656  this->envConst_.getRawPtr()
657  //,task_communication_xadj,
658  //task_communication_adj
659  );
660 
661  // KDD For now, we would need to re-map the part numbers in the solution.
662  // KDD I suspect we'll later need to distinguish between part numbers and
663  // KDD process numbers to provide separation between partitioning and
664  // KDD mapping. For example, does this approach here assume #parts == #procs?
665  // KDD If we map k tasks to p processes with k > p, do we effectively reduce
666  // KDD the number of tasks (parts) in the solution?
667 
668 #ifdef KDD_READY
669  const part_t *oldParts = solution_->getPartListView();
670  size_t nLocal = ia->getNumLocalIds();
671  for (size_t i = 0; i < nLocal; i++) {
672  // kind of cheating since oldParts is a view; probably want an interface in solution
673  // for resetting the PartList rather than hacking in like this.
674  oldParts[i] = ctm->getAssignedProcForTask(oldParts[i]);
675  }
676 #endif
677 
678  //for now just delete the object.
679  delete ctm;
680  }
681 #endif
682 
683  else if (mapping_type == 1){
684  //if mapping is 1 -- graph mapping
685  }
686 
687  this->env_->debug(DETAILED_STATUS, "Exiting solve");
688 }
689 
690 template <typename Adapter>
692 {
693  this->env_->debug(DETAILED_STATUS,
694  "PartitioningProblem::createPartitioningProblem");
695 
696  using Teuchos::ParameterList;
697 
698  // A Problem object may be reused. The input data may have changed and
699  // new parameters or part sizes may have been set.
700  //
701  // Save these values in order to determine if we need to create a new model.
702 
703  //ModelType previousModel = modelType_;
704  bool prevModelAvail[MAX_NUM_MODEL_TYPES];
705  for(int i=0;i<MAX_NUM_MODEL_TYPES;i++)
706  {
707  prevModelAvail[i] = modelAvail_[i];
708  }
709 
710 
711  modelFlag_t previousGraphModelFlags = graphFlags_;
712  modelFlag_t previousIdentifierModelFlags = idFlags_;
713  modelFlag_t previousCoordinateModelFlags = coordFlags_;
714 
715  //modelType_ = InvalidModel;
716  for(int i=0;i<MAX_NUM_MODEL_TYPES;i++)
717  {
718  modelAvail_[i] = false;
719  }
720 
721  graphFlags_.reset();
722  idFlags_.reset();
723  coordFlags_.reset();
724 
726  // It's possible at this point that the Problem may want to
727  // add problem parameters to the parameter list in the Environment.
728  //
729  // Since the parameters in the Environment have already been
730  // validated in its constructor, a new Environment must be created:
732  // Teuchos::RCP<const Teuchos::Comm<int> > oldComm = this->env_->comm_;
733  // const ParameterList &oldParams = this->env_->getUnvalidatedParameters();
734  //
735  // ParameterList newParams = oldParams;
736  // newParams.set("new_parameter", "new_value");
737  //
738  // ParameterList &newPartParams = newParams.sublist("partitioning");
739  // newPartParams.set("new_partitioning_parameter", "its_value");
740  //
741  // this->env_ = rcp(new Environment(newParams, oldComm));
743 
744  this->env_->debug(DETAILED_STATUS, " parameters");
745  Environment &env = *(this->env_);
746  ParameterList &pl = env.getParametersNonConst();
747 
748  std::string defString("default");
749 
750  // Did the user specify a computational model?
751 
752  std::string model(defString);
753  const Teuchos::ParameterEntry *pe = pl.getEntryPtr("model");
754  if (pe)
755  model = pe->getValue<std::string>(&model);
756 
757  // Did the user specify an algorithm?
758 
759  std::string algorithm(defString);
760  pe = pl.getEntryPtr("algorithm");
761  if (pe)
762  algorithm = pe->getValue<std::string>(&algorithm);
763 
764  // Possible algorithm requirements that must be conveyed to the model:
765 
766  bool needConsecutiveGlobalIds = false;
767  bool removeSelfEdges= false;
768 
770  // Determine algorithm, model, and algorithm requirements. This
771  // is a first pass. Feel free to change this and add to it.
772 
773  if (algorithm != defString)
774  {
775 
776  // Figure out the model required by the algorithm
777  if (algorithm == std::string("block") ||
778  algorithm == std::string("random") ||
779  algorithm == std::string("cyclic") ){
780 
781  //modelType_ = IdentifierModelType;
782  modelAvail_[IdentifierModelType] = true;
783 
784  algName_ = algorithm;
785  }
786  else if (algorithm == std::string("zoltan") ||
787  algorithm == std::string("parma") ||
788  algorithm == std::string("forTestingOnly"))
789  {
790  algName_ = algorithm;
791  }
792  else if (algorithm == std::string("rcb") ||
793  algorithm == std::string("rib") ||
794  algorithm == std::string("hsfc"))
795  {
796  // rcb, rib, hsfc provided through Zoltan
797  Teuchos::ParameterList &zparams = pl.sublist("zoltan_parameters",false);
798  zparams.set("LB_METHOD", algorithm);
799  if (numberOfWeights_ > 0) {
800  char strval[20];
801  sprintf(strval, "%d", numberOfWeights_);
802  zparams.set("OBJ_WEIGHT_DIM", strval);
803  }
804  algName_ = std::string("zoltan");
805  }
806  else if (algorithm == std::string("multijagged"))
807  {
808  //modelType_ = CoordinateModelType;
809  modelAvail_[CoordinateModelType]=true;
810 
811  algName_ = algorithm;
812  }
813  else if (algorithm == std::string("metis") ||
814  algorithm == std::string("parmetis"))
815  {
816 
817  //modelType_ = GraphModelType;
818  modelAvail_[GraphModelType]=true;
819  algName_ = algorithm;
820  removeSelfEdges = true;
821  needConsecutiveGlobalIds = true;
822  }
823  else if (algorithm == std::string("scotch") ||
824  algorithm == std::string("ptscotch")) // BDD: Don't construct graph for scotch here
825  {
826  algName_ = algorithm;
827  }
828  else if (algorithm == std::string("pulp"))
829  {
830  algName_ = algorithm;
831  }
832  else if (algorithm == std::string("patoh") ||
833  algorithm == std::string("phg"))
834  {
835  // if ((modelType_ != GraphModelType) &&
836  // (modelType_ != HypergraphModelType) )
837  if ((modelAvail_[GraphModelType]==false) &&
838  (modelAvail_[HypergraphModelType]==false) )
839  {
840  //modelType_ = HypergraphModelType;
841  modelAvail_[HypergraphModelType]=true;
842  }
843  algName_ = algorithm;
844  }
845  else
846  {
847  // Parameter list should ensure this does not happen.
848  throw std::logic_error("parameter list algorithm is invalid");
849  }
850  }
851  else if (model != defString)
852  {
853  // Figure out the algorithm suggested by the model.
854  if (model == std::string("hypergraph"))
855  {
856  //modelType_ = HypergraphModelType;
857  modelAvail_[HypergraphModelType]=true;
858 
859  algName_ = std::string("phg");
860  }
861  else if (model == std::string("graph"))
862  {
863  //modelType_ = GraphModelType;
864  modelAvail_[GraphModelType]=true;
865 
866 #ifdef HAVE_ZOLTAN2_SCOTCH
867  modelAvail_[GraphModelType]=false; // graph constructed by AlgPTScotch
868  if (this->comm_->getSize() > 1)
869  algName_ = std::string("ptscotch");
870  else
871  algName_ = std::string("scotch");
872 #else
873 #ifdef HAVE_ZOLTAN2_PARMETIS
874  if (this->comm_->getSize() > 1)
875  algName_ = std::string("parmetis");
876  else
877  algName_ = std::string("metis");
878  removeSelfEdges = true;
879  needConsecutiveGlobalIds = true;
880 #else
881 #ifdef HAVE_ZOLTAN2_PULP
882  // TODO: XtraPuLP
883  //if (this->comm_->getSize() > 1)
884  // algName_ = std::string("xtrapulp");
885  //else
886  algName_ = std::string("pulp");
887 #else
888  algName_ = std::string("phg");
889 #endif
890 #endif
891 #endif
892  }
893  else if (model == std::string("geometry"))
894  {
895  //modelType_ = CoordinateModelType;
896  modelAvail_[CoordinateModelType]=true;
897 
898  algName_ = std::string("multijagged");
899  }
900  else if (model == std::string("ids"))
901  {
902  //modelType_ = IdentifierModelType;
903  modelAvail_[IdentifierModelType]=true;
904 
905  algName_ = std::string("block");
906  }
907  else
908  {
909  // Parameter list should ensure this does not happen.
910  env.localBugAssertion(__FILE__, __LINE__,
911  "parameter list model type is invalid", 1, BASIC_ASSERTION);
912  }
913  }
914  else
915  {
916  // Determine an algorithm and model suggested by the input type.
917  // TODO: this is a good time to use the time vs. quality parameter
918  // in choosing an algorithm, and setting some parameters
919 
920  if (inputType_ == MatrixAdapterType)
921  {
922  //modelType_ = HypergraphModelType;
923  modelAvail_[HypergraphModelType]=true;
924 
925  algName_ = std::string("phg");
926  }
927  else if (inputType_ == GraphAdapterType ||
928  inputType_ == MeshAdapterType)
929  {
930  //modelType_ = GraphModelType;
931  modelAvail_[GraphModelType]=true;
932 
933  algName_ = std::string("phg");
934  }
935  else if (inputType_ == VectorAdapterType)
936  {
937  //modelType_ = CoordinateModelType;
938  modelAvail_[CoordinateModelType]=true;
939 
940  algName_ = std::string("multijagged");
941  }
942  else if (inputType_ == IdentifierAdapterType)
943  {
944  //modelType_ = IdentifierModelType;
945  modelAvail_[IdentifierModelType]=true;
946 
947  algName_ = std::string("block");
948  }
949  else{
950  // This should never happen
951  throw std::logic_error("input type is invalid");
952  }
953  }
954 
955  // Hierarchical partitioning?
956 
957  Array<int> valueList;
958  pe = pl.getEntryPtr("topology");
959 
960  if (pe){
961  valueList = pe->getValue<Array<int> >(&valueList);
962 
963  if (!Zoltan2::noValuesAreInRangeList<int>(valueList)){
964  int *n = new int [valueList.size() + 1];
965  levelNumberParts_ = arcp(n, 0, valueList.size() + 1, true);
966  int procsPerNode = 1;
967  for (int i=0; i < valueList.size(); i++){
968  levelNumberParts_[i+1] = valueList[i];
969  procsPerNode *= valueList[i];
970  }
971  // Number of parts in the first level
972  levelNumberParts_[0] = env.numProcs_ / procsPerNode;
973 
974  if (env.numProcs_ % procsPerNode > 0)
975  levelNumberParts_[0]++;
976  }
977  }
978  else{
979  levelNumberParts_.clear();
980  }
981 
982  hierarchical_ = levelNumberParts_.size() > 0;
983 
984  // Object to be partitioned? (rows, columns, etc)
985 
986  std::string objectOfInterest(defString);
987  pe = pl.getEntryPtr("objects_to_partition");
988  if (pe)
989  objectOfInterest = pe->getValue<std::string>(&objectOfInterest);
990 
992  // Set model creation flags, if any.
993 
994  this->env_->debug(DETAILED_STATUS, " models");
995  // if (modelType_ == GraphModelType)
996  if (modelAvail_[GraphModelType]==true)
997  {
998 
999  // Any parameters in the graph sublist?
1000 
1001  std::string symParameter(defString);
1002  pe = pl.getEntryPtr("symmetrize_graph");
1003  if (pe){
1004  symParameter = pe->getValue<std::string>(&symParameter);
1005  if (symParameter == std::string("transpose"))
1006  graphFlags_.set(SYMMETRIZE_INPUT_TRANSPOSE);
1007  else if (symParameter == std::string("bipartite"))
1008  graphFlags_.set(SYMMETRIZE_INPUT_BIPARTITE);
1009  }
1010 
1011  bool sgParameter = false;
1012  pe = pl.getEntryPtr("subset_graph");
1013  if (pe)
1014  sgParameter = pe->getValue(&sgParameter);
1015 
1016  if (sgParameter == 1)
1017  graphFlags_.set(BUILD_SUBSET_GRAPH);
1018 
1019  // Any special behaviors required by the algorithm?
1020 
1021  if (removeSelfEdges)
1022  graphFlags_.set(REMOVE_SELF_EDGES);
1023 
1024  if (needConsecutiveGlobalIds)
1025  graphFlags_.set(GENERATE_CONSECUTIVE_IDS);
1026 
1027  // How does user input map to vertices and edges?
1028 
1029  if (inputType_ == MatrixAdapterType){
1030  if (objectOfInterest == defString ||
1031  objectOfInterest == std::string("matrix_rows") )
1032  graphFlags_.set(VERTICES_ARE_MATRIX_ROWS);
1033  else if (objectOfInterest == std::string("matrix_columns"))
1034  graphFlags_.set(VERTICES_ARE_MATRIX_COLUMNS);
1035  else if (objectOfInterest == std::string("matrix_nonzeros"))
1036  graphFlags_.set(VERTICES_ARE_MATRIX_NONZEROS);
1037  }
1038 
1039  else if (inputType_ == MeshAdapterType){
1040  if (objectOfInterest == defString ||
1041  objectOfInterest == std::string("mesh_nodes") )
1042  graphFlags_.set(VERTICES_ARE_MESH_NODES);
1043  else if (objectOfInterest == std::string("mesh_elements"))
1044  graphFlags_.set(VERTICES_ARE_MESH_ELEMENTS);
1045  }
1046  }
1047  //MMW is it ok to remove else?
1048  // else if (modelType_ == IdentifierModelType)
1049  if (modelAvail_[IdentifierModelType]==true)
1050  {
1051 
1052  // Any special behaviors required by the algorithm?
1053 
1054  }
1055  // else if (modelType_ == CoordinateModelType)
1056  if (modelAvail_[CoordinateModelType]==true)
1057  {
1058 
1059  // Any special behaviors required by the algorithm?
1060 
1061  }
1062 
1063 
1064  if (newData ||
1065  (modelAvail_[GraphModelType]!=prevModelAvail[GraphModelType]) ||
1066  (modelAvail_[HypergraphModelType]!=prevModelAvail[HypergraphModelType])||
1067  (modelAvail_[CoordinateModelType]!=prevModelAvail[CoordinateModelType])||
1068  (modelAvail_[IdentifierModelType]!=prevModelAvail[IdentifierModelType])||
1069  // (modelType_ != previousModel) ||
1070  (graphFlags_ != previousGraphModelFlags) ||
1071  (coordFlags_ != previousCoordinateModelFlags) ||
1072  (idFlags_ != previousIdentifierModelFlags))
1073  {
1074  // Create the computational model.
1075  // Models are instantiated for base input adapter types (mesh,
1076  // matrix, graph, and so on). We pass a pointer to the input
1077  // adapter, cast as the base input type.
1078 
1079  //KDD Not sure why this shadow declaration is needed
1080  //KDD Comment out for now; revisit later if problems.
1081  //KDD const Teuchos::ParameterList pl = this->envConst_->getParameters();
1082  //bool exceptionThrow = true;
1083 
1084  if(modelAvail_[GraphModelType]==true)
1085  {
1086  this->env_->debug(DETAILED_STATUS, " building graph model");
1087  this->graphModel_ = rcp(new GraphModel<base_adapter_t>(
1088  this->baseInputAdapter_, this->envConst_, this->comm_,
1089  graphFlags_));
1090 
1091  this->baseModel_ = rcp_implicit_cast<const Model<base_adapter_t> >(
1092  this->graphModel_);
1093  }
1094  if(modelAvail_[HypergraphModelType]==true)
1095  {
1096  //KDD USING ZOLTAN FOR HYPERGRAPH FOR NOW
1097  //KDD std::cout << "Hypergraph model not implemented yet..." << std::endl;
1098  }
1099 
1100  if(modelAvail_[CoordinateModelType]==true)
1101  {
1102  this->env_->debug(DETAILED_STATUS, " building coordinate model");
1103  this->coordinateModel_ = rcp(new CoordinateModel<base_adapter_t>(
1104  this->baseInputAdapter_, this->envConst_, this->comm_,
1105  coordFlags_));
1106 
1107  this->baseModel_ = rcp_implicit_cast<const Model<base_adapter_t> >(
1108  this->coordinateModel_);
1109  }
1110 
1111  if(modelAvail_[IdentifierModelType]==true)
1112  {
1113  this->env_->debug(DETAILED_STATUS, " building identifier model");
1114  this->identifierModel_ = rcp(new IdentifierModel<base_adapter_t>(
1115  this->baseInputAdapter_, this->envConst_, this->comm_,
1116  idFlags_));
1117 
1118  this->baseModel_ = rcp_implicit_cast<const Model<base_adapter_t> >(
1119  this->identifierModel_);
1120  }
1121 
1122  this->env_->memory("After creating Model");
1123  this->env_->debug(DETAILED_STATUS, "createPartitioningProblem done");
1124  }
1125 }
1126 
1127 } // namespace Zoltan2
1128 #endif
use columns as graph vertices
Serial greedy first-fit graph coloring (local graph only)
#define HELLO
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Time an algorithm (or other entity) as a whole.
ignore invalid neighbors
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
PartitioningProblem(Adapter *A, ParameterList *p, const RCP< const Teuchos::Comm< int > > &comm)
Constructor where Teuchos communicator is specified.
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&#39;s entry and exit
static void getValidParameters(ParameterList &pl)
Set up validators specific to this Problem.
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.
SparseMatrixAdapter_t::part_t part_t
use matrix rows as graph vertices
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.
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.
Defines the Problem base class.
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.
Define IntegerRangeList validator.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
use mesh elements as vertices
Defines the GraphModel interface.
static void getValidParameters(ParameterList &pl)
Set up validators specific to this algorithm.
static void getValidParameters(ParameterList &)
Set up validators specific to this algorithm.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Multi Jagged coordinate partitioning algorithm.
static void getValidParameters(ParameterList &pl)
Set up validators specific to this algorithm.
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:74