Zoltan2
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>
61 #include <Zoltan2_TaskMapping.hpp>
62 
63 #ifndef _WIN32
64 #include <unistd.h>
65 #else
66 #include <process.h>
67 #define NOMINMAX
68 #include <windows.h>
69 #endif
70 
71 #ifdef HAVE_ZOLTAN2_OVIS
72 #include <ovis.h>
73 #endif
74 
75 namespace Zoltan2{
76 
100 template<typename Adapter>
101 class PartitioningProblem : public Problem<Adapter>
102 {
103 public:
104 
105  typedef typename Adapter::scalar_t scalar_t;
106  typedef typename Adapter::gno_t gno_t;
107  typedef typename Adapter::lno_t lno_t;
108  typedef typename Adapter::part_t part_t;
109  typedef typename Adapter::user_t user_t;
110  typedef typename Adapter::base_adapter_t base_adapter_t;
111 
112 #ifdef HAVE_ZOLTAN2_MPI
113  typedef Teuchos::OpaqueWrapper<MPI_Comm> mpiWrapper_t;
116  PartitioningProblem(Adapter *A, ParameterList *p, MPI_Comm comm):
117  Problem<Adapter>(A,p,comm), solution_(),
118  problemComm_(), problemCommConst_(),
119  inputType_(InvalidAdapterType),
120  graphFlags_(), idFlags_(), coordFlags_(), algName_(),
121  numberOfWeights_(), partIds_(), partSizes_(),
122  numberOfCriteria_(), levelNumberParts_(), hierarchical_(false),
123  metricsRequested_(false), graphMetricsRequested_(false), metrics_(),
124  graphMetrics_()
125  {
126  for(int i=0;i<MAX_NUM_MODEL_TYPES;i++) modelAvail_[i]=false;
127  initializeProblem();
128  }
129 #endif
130 
132  PartitioningProblem(Adapter *A, ParameterList *p):
133  Problem<Adapter>(A,p), solution_(),
134  problemComm_(), problemCommConst_(),
135  inputType_(InvalidAdapterType),
136  graphFlags_(), idFlags_(), coordFlags_(), algName_(),
137  numberOfWeights_(),
138  partIds_(), partSizes_(), numberOfCriteria_(),
139  levelNumberParts_(), hierarchical_(false),
140  metricsRequested_(false), graphMetricsRequested_(false), metrics_(),
141  graphMetrics_()
142  {
143  for(int i=0;i<MAX_NUM_MODEL_TYPES;i++) modelAvail_[i]=false;
144  initializeProblem();
145  }
146 
148  PartitioningProblem(Adapter *A, ParameterList *p,
149  RCP<const Teuchos::Comm<int> > &comm):
150  Problem<Adapter>(A,p,comm), solution_(),
151  problemComm_(), problemCommConst_(),
152  inputType_(InvalidAdapterType),
153  graphFlags_(), idFlags_(), coordFlags_(), algName_(),
154  numberOfWeights_(),
155  partIds_(), partSizes_(), numberOfCriteria_(),
156  levelNumberParts_(), hierarchical_(false),
157  metricsRequested_(false), graphMetricsRequested_(false), metrics_(),
158  graphMetrics_()
159  {
160  for(int i=0;i<MAX_NUM_MODEL_TYPES;i++) modelAvail_[i]=false;
161  initializeProblem();
162  }
163 
167 
169  //
170  // \param updateInputData If true this indicates that either
171  // this is the first attempt at solution, or that we
172  // are computing a new solution and the input data has
173  // changed since the previous solution was computed.
174  // By input data we mean coordinates, topology, or weights.
175  // If false, this indicates that we are computing a
176  // new solution using the same input data was used for
177  // the previous solution, even though the parameters
178  // may have been changed.
179  //
180  // For the sake of performance, we ask the caller to set \c updateInputData
181  // to false if he/she is computing a new solution using the same input data,
182  // but different problem parameters, than that which was used to compute
183  // the most recent solution.
184 
185  void solve(bool updateInputData=true );
186 
188  //
189  // \return a reference to the solution to the most recent solve().
190 
192  return *(solution_.getRawPtr());
193  };
194 
203  const scalar_t getWeightImbalance(int idx=0) const {
204  scalar_t imb = 0;
205  if (!metrics_.is_null())
206  metrics_->getWeightImbalance(imb, idx);
207 
208  return imb;
209  }
210 
215  ArrayRCP<const MetricValues<scalar_t> > getMetrics() const {
216  if (metrics_.is_null()){
217  ArrayRCP<const MetricValues<scalar_t> > emptyMetrics;
218  return emptyMetrics;
219  }
220  else
221  return metrics_->getMetrics();
222  }
223 
228  ArrayRCP<const GraphMetricValues<scalar_t> > getGraphMetrics() const {
229  if (graphMetrics_.is_null()){
230  ArrayRCP<const GraphMetricValues<scalar_t> > emptyMetrics;
231  return emptyMetrics;
232  }
233  else
234  return graphMetrics_->getGraphMetrics();
235  }
236 
242  void printMetrics(std::ostream &os) const {
243  if (metrics_.is_null())
244  os << "No metrics available." << std::endl;
245  else
246  metrics_->printMetrics(os);
247  };
248 
287  void setPartSizes(int len, part_t *partIds, scalar_t *partSizes,
288  bool makeCopy=true)
289  {
290  setPartSizesForCriteria(0, len, partIds, partSizes, makeCopy);
291  }
292 
327  void setPartSizesForCriteria(int criteria, int len, part_t *partIds,
328  scalar_t *partSizes, bool makeCopy=true) ;
329 /*
330  void setMachine(MachineRepresentation<typename Adapter::base_adapter_t::scalar_t> *machine);
331 */
334  void resetParameters(ParameterList *params)
335  {
336  Problem<Adapter>::resetParameters(params); // creates new environment
337  }
338 
343  const RCP<const Environment> & getEnvironment() const
344  {
345  return this->envConst_;
346  }
347 
348 private:
349  void initializeProblem();
350 
351  void createPartitioningProblem(bool newData);
352 
353  RCP<PartitioningSolution<Adapter> > solution_;
354 
355  RCP<MachineRepresentation <typename Adapter::base_adapter_t::scalar_t> > machine_;
356 
357  RCP<Comm<int> > problemComm_;
358  RCP<const Comm<int> > problemCommConst_;
359 
360  BaseAdapterType inputType_;
361 
362  //ModelType modelType_;
363  bool modelAvail_[MAX_NUM_MODEL_TYPES];
364 
365  modelFlag_t graphFlags_;
366  modelFlag_t idFlags_;
367  modelFlag_t coordFlags_;
368  std::string algName_;
369 
370  int numberOfWeights_;
371 
372  // Suppose Array<part_t> partIds = partIds_[w]. If partIds.size() > 0
373  // then the user supplied part sizes for weight index "w", and the sizes
374  // corresponding to the Ids in partIds are partSizes[w].
375  //
376  // If numberOfWeights_ >= 0, then there is an Id and Sizes array for
377  // for each weight. Otherwise the user did not supply object weights,
378  // but they can still specify part sizes.
379  // So numberOfCriteria_ is numberOfWeights_ or one, whichever is greater.
380 
381  ArrayRCP<ArrayRCP<part_t> > partIds_;
382  ArrayRCP<ArrayRCP<scalar_t> > partSizes_;
383  int numberOfCriteria_;
384 
385  // Number of parts to be computed at each level in hierarchical partitioning.
386 
387  ArrayRCP<int> levelNumberParts_;
388  bool hierarchical_;
389 
390  // Did the user request metrics?
391 
392  bool metricsRequested_;
393  bool graphMetricsRequested_;
394  RCP<const PartitioningSolutionQuality<Adapter> > metrics_;
395  RCP<const GraphPartitioningSolutionQuality<Adapter> > graphMetrics_;
396 };
398 
399 /*
400 template <typename Adapter>
401 void PartitioningProblem<Adapter>::setMachine(MachineRepresentation<typename Adapter::base_adapter_t::scalar_t> *machine){
402  this->machine_ = RCP<MachineRepresentation<typename Adapter::base_adapter_t::scalar_t> > (machine, false);
403 }
404 */
405 
406 
407 template <typename Adapter>
409 {
410  HELLO;
411 
412  this->env_->debug(DETAILED_STATUS, "PartitioningProblem::initializeProblem");
413 
414  if (getenv("DEBUGME")){
415 #ifndef _WIN32
416  std::cout << getpid() << std::endl;
417  sleep(15);
418 #else
419  std::cout << _getpid() << std::endl;
420  Sleep(15000);
421 #endif
422  }
423 
424 #ifdef HAVE_ZOLTAN2_OVIS
425  ovis_enabled(this->comm_->getRank());
426 #endif
427 
428  // Create a copy of the user's communicator.
429 
430  problemComm_ = this->comm_->duplicate();
431  problemCommConst_ = rcp_const_cast<const Comm<int> > (problemComm_);
432 
433  machine_ = RCP <Zoltan2::MachineRepresentation<typename Adapter::scalar_t> >(new Zoltan2::MachineRepresentation<typename Adapter::scalar_t>(problemComm_));
434 
435  // Number of criteria is number of user supplied weights if non-zero.
436  // Otherwise it is 1 and uniform weight is implied.
437 
438  numberOfWeights_ = this->inputAdapter_->getNumWeightsPerID();
439 
440  numberOfCriteria_ = (numberOfWeights_ > 1) ? numberOfWeights_ : 1;
441 
442  inputType_ = this->inputAdapter_->adapterType();
443 
444  // The Caller can specify part sizes in setPartSizes(). If he/she
445  // does not, the part size arrays are empty.
446 
447  ArrayRCP<part_t> *noIds = new ArrayRCP<part_t> [numberOfCriteria_];
448  ArrayRCP<scalar_t> *noSizes = new ArrayRCP<scalar_t> [numberOfCriteria_];
449 
450  partIds_ = arcp(noIds, 0, numberOfCriteria_, true);
451  partSizes_ = arcp(noSizes, 0, numberOfCriteria_, true);
452 
453  if (this->env_->getDebugLevel() >= DETAILED_STATUS){
454  std::ostringstream msg;
455  msg << problemComm_->getSize() << " procs,"
456  << numberOfWeights_ << " user-defined weights\n";
457  this->env_->debug(DETAILED_STATUS, msg.str());
458  }
459 
460  this->env_->memory("After initializeProblem");
461 }
462 
463 template <typename Adapter>
465  int criteria, int len, part_t *partIds, scalar_t *partSizes, bool makeCopy)
466 {
467  this->env_->localInputAssertion(__FILE__, __LINE__, "invalid length",
468  len>= 0, BASIC_ASSERTION);
469 
470  this->env_->localInputAssertion(__FILE__, __LINE__, "invalid criteria",
471  criteria >= 0 && criteria < numberOfCriteria_, BASIC_ASSERTION);
472 
473  if (len == 0){
474  partIds_[criteria] = ArrayRCP<part_t>();
475  partSizes_[criteria] = ArrayRCP<scalar_t>();
476  return;
477  }
478 
479  this->env_->localInputAssertion(__FILE__, __LINE__, "invalid arrays",
480  partIds && partSizes, BASIC_ASSERTION);
481 
482  // The global validity of the partIds and partSizes arrays is performed
483  // by the PartitioningSolution, which computes global part distribution and
484  // part sizes.
485 
486  part_t *z2_partIds = NULL;
487  scalar_t *z2_partSizes = NULL;
488  bool own_memory = false;
489 
490  if (makeCopy){
491  z2_partIds = new part_t [len];
492  z2_partSizes = new scalar_t [len];
493  this->env_->localMemoryAssertion(__FILE__, __LINE__, len, z2_partSizes);
494  memcpy(z2_partIds, partIds, len * sizeof(part_t));
495  memcpy(z2_partSizes, partSizes, len * sizeof(scalar_t));
496  own_memory=true;
497  }
498  else{
499  z2_partIds = partIds;
500  z2_partSizes = partSizes;
501  }
502 
503  partIds_[criteria] = arcp(z2_partIds, 0, len, own_memory);
504  partSizes_[criteria] = arcp(z2_partSizes, 0, len, own_memory);
505 }
506 
507 template <typename Adapter>
508 void PartitioningProblem<Adapter>::solve(bool updateInputData)
509 {
510  this->env_->debug(DETAILED_STATUS, "Entering solve");
511 
512  // Create the computational model.
513 
514  this->env_->timerStart(MACRO_TIMERS, "create problem");
515 
516  createPartitioningProblem(updateInputData);
517 
518  this->env_->timerStop(MACRO_TIMERS, "create problem");
519 
520  // TODO: If hierarchical_
521 
522  // Create the solution. The algorithm will query the Solution
523  // for part and weight information. The algorithm will
524  // update the solution with part assignments and quality metrics.
525 
526  // Create the algorithm
527  try {
528  if (algName_ == std::string("multijagged")) {
529  this->algorithm_ = rcp(new Zoltan2_AlgMJ<Adapter>(this->envConst_,
530  problemComm_,
531  this->coordinateModel_));
532  }
533  else if (algName_ == std::string("zoltan")) {
534  this->algorithm_ = rcp(new AlgZoltan<Adapter>(this->envConst_,
535  problemComm_,
536  this->baseInputAdapter_));
537  }
538  else if (algName_ == std::string("parma")) {
539  this->algorithm_ = rcp(new AlgParMA<Adapter>(this->envConst_,
540  problemComm_,
541  this->baseInputAdapter_));
542  }
543  else if (algName_ == std::string("scotch")) {
544  this->algorithm_ = rcp(new AlgPTScotch<Adapter>(this->envConst_,
545  problemComm_,
546  this->graphModel_));
547  }
548  else if (algName_ == std::string("parmetis")) {
549  this->algorithm_ = rcp(new AlgParMETIS<Adapter>(this->envConst_,
550  problemComm_,
551  this->graphModel_));
552  }
553  else if (algName_ == std::string("block")) {
554  this->algorithm_ = rcp(new AlgBlock<Adapter>(this->envConst_,
555  problemComm_, this->identifierModel_));
556  }
557  else if (algName_ == std::string("forTestingOnly")) {
558  this->algorithm_ = rcp(new AlgForTestingOnly<Adapter>(this->envConst_,
559  problemComm_,
560  this->baseInputAdapter_));
561  }
562  // else if (algName_ == std::string("rcb")) {
563  // this->algorithm_ = rcp(new AlgRCB<Adapter>(this->envConst_,problemComm_,
564  // this->coordinateModel_));
565  // }
566  else {
567  throw std::logic_error("partitioning algorithm not supported");
568  }
569  }
571 
572  // Create the solution
573  this->env_->timerStart(MACRO_TIMERS, "create solution");
574  PartitioningSolution<Adapter> *soln = NULL;
575 
576  try{
577  soln = new PartitioningSolution<Adapter>(
578  this->envConst_, problemCommConst_, numberOfWeights_,
579  partIds_.view(0, numberOfCriteria_),
580  partSizes_.view(0, numberOfCriteria_), this->algorithm_);
581  }
583 
584  solution_ = rcp(soln);
585 
586  this->env_->timerStop(MACRO_TIMERS, "create solution");
587  this->env_->memory("After creating Solution");
588 
589  // Call the algorithm
590 
591  try {
592  this->algorithm_->partition(solution_);
593  }
595 
596  //if mapping is requested
597  const Teuchos::ParameterEntry *pe = this->envConst_->getParameters().getEntryPtr("mapping_type");
598  int mapping_type = -1;
599  if (pe){
600  mapping_type = pe->getValue(&mapping_type);
601  }
602  //if mapping is 0 -- coordinate mapping
603 
604  if (mapping_type == 0){
605 
606  //part_t *task_communication_xadj = NULL, *task_communication_adj = NULL;
607 
610  problemComm_.getRawPtr(),
611  machine_.getRawPtr(),
612  this->coordinateModel_.getRawPtr(),
613  solution_.getRawPtr(),
614  this->envConst_.getRawPtr()
615  //,task_communication_xadj,
616  //task_communication_adj
617  );
618  //for now just delete the object.
619  delete ctm;
620  }
621  else if (mapping_type == 1){
622  //if mapping is 1 -- graph mapping
623  }
624 
625  if (metricsRequested_){
626  typedef PartitioningSolution<Adapter> ps_t;
628 
629  psq_t *quality = NULL;
630  RCP<const ps_t> solutionConst = rcp_const_cast<const ps_t>(solution_);
631 
632  try{
633  quality = new psq_t(this->envConst_, problemCommConst_,
634  this->inputAdapter_, solutionConst);
635  }
636  Z2_FORWARD_EXCEPTIONS
637 
638  metrics_ = rcp(quality);
639  }
640 
641  if (graphMetricsRequested_ && (algName_ == std::string("scotch") ||
642  algName_ == std::string("parmetis"))){
643  typedef PartitioningSolution<Adapter> ps_t;
645 
646  gpsq_t *quality = NULL;
647  RCP<const ps_t> solutionConst = rcp_const_cast<const ps_t>(solution_);
648 
649  try{
650  quality = new gpsq_t(this->envConst_, problemCommConst_,
651  this->graphModel_, this->inputAdapter_,
652  solutionConst);
653  }
654  Z2_FORWARD_EXCEPTIONS
655 
656  graphMetrics_ = rcp(quality);
657  }
658 
659  this->env_->debug(DETAILED_STATUS, "Exiting solve");
660 }
661 
662 template <typename Adapter>
664 {
665  this->env_->debug(DETAILED_STATUS,
666  "PartitioningProblem::createPartitioningProblem");
667 
668  using Teuchos::ParameterList;
669 
670  // A Problem object may be reused. The input data may have changed and
671  // new parameters or part sizes may have been set.
672  //
673  // Save these values in order to determine if we need to create a new model.
674 
675  //ModelType previousModel = modelType_;
676  bool prevModelAvail[MAX_NUM_MODEL_TYPES];
677  for(int i=0;i<MAX_NUM_MODEL_TYPES;i++)
678  {
679  prevModelAvail[i] = modelAvail_[i];
680  }
681 
682 
683  modelFlag_t previousGraphModelFlags = graphFlags_;
684  modelFlag_t previousIdentifierModelFlags = idFlags_;
685  modelFlag_t previousCoordinateModelFlags = coordFlags_;
686 
687  //modelType_ = InvalidModel;
688  for(int i=0;i<MAX_NUM_MODEL_TYPES;i++)
689  {
690  modelAvail_[i] = false;
691  }
692 
693  graphFlags_.reset();
694  idFlags_.reset();
695  coordFlags_.reset();
696 
698  // It's possible at this point that the Problem may want to
699  // add problem parameters to the parameter list in the Environment.
700  //
701  // Since the parameters in the Environment have already been
702  // validated in its constructor, a new Environment must be created:
704  // Teuchos::RCP<const Teuchos::Comm<int> > oldComm = this->env_->comm_;
705  // const ParameterList &oldParams = this->env_->getUnvalidatedParameters();
706  //
707  // ParameterList newParams = oldParams;
708  // newParams.set("new_parameter", "new_value");
709  //
710  // ParameterList &newPartParams = newParams.sublist("partitioning");
711  // newPartParams.set("new_partitioning_parameter", "its_value");
712  //
713  // this->env_ = rcp(new Environment(newParams, oldComm));
715 
716  this->env_->debug(DETAILED_STATUS, " parameters");
717  Environment &env = *(this->env_);
718  ParameterList &pl = env.getParametersNonConst();
719 
720  std::string defString("default");
721 
722  // Did the user ask for computation of quality metrics?
723 
724  int yesNo=0;
725  const Teuchos::ParameterEntry *pe = pl.getEntryPtr("compute_metrics");
726  if (pe){
727  yesNo = pe->getValue<int>(&yesNo);
728  metricsRequested_ = true;
729  }
730 
731  // Did the user specify a computational model?
732 
733  std::string model(defString);
734  pe = pl.getEntryPtr("model");
735  if (pe)
736  model = pe->getValue<std::string>(&model);
737 
738  // Did the user specify an algorithm?
739 
740  std::string algorithm(defString);
741  pe = pl.getEntryPtr("algorithm");
742  if (pe)
743  algorithm = pe->getValue<std::string>(&algorithm);
744 
745  // Possible algorithm requirements that must be conveyed to the model:
746 
747  bool needConsecutiveGlobalIds = false;
748  bool removeSelfEdges= false;
749 
751  // Determine algorithm, model, and algorithm requirements. This
752  // is a first pass. Feel free to change this and add to it.
753 
754  if (algorithm != defString)
755  {
756 
757  // Figure out the model required by the algorithm
758  if (algorithm == std::string("block") ||
759  algorithm == std::string("random") ||
760  algorithm == std::string("cyclic") ){
761 
762  //modelType_ = IdentifierModelType;
763  modelAvail_[IdentifierModelType] = true;
764 
765  algName_ = algorithm;
766  }
767  else if (algorithm == std::string("zoltan") ||
768  algorithm == std::string("parma") ||
769  algorithm == std::string("forTestingOnly"))
770  {
771  algName_ = algorithm;
772  }
773  else if (algorithm == std::string("rcb") ||
774  algorithm == std::string("rib") ||
775  algorithm == std::string("hsfc"))
776  {
777  // rcb, rib, hsfc provided through Zoltan
778  Teuchos::ParameterList &zparams = pl.sublist("zoltan_parameters",false);
779  zparams.set("LB_METHOD", algorithm);
780  if (numberOfWeights_ > 0) {
781  char strval[10];
782  sprintf(strval, "%d", numberOfWeights_);
783  zparams.set("OBJ_WEIGHT_DIM", strval);
784  }
785  algName_ = std::string("zoltan");
786  }
787  else if (algorithm == std::string("multijagged"))
788  {
789  //modelType_ = CoordinateModelType;
790  modelAvail_[CoordinateModelType]=true;
791 
792  algName_ = algorithm;
793  }
794  else if (algorithm == std::string("metis") ||
795  algorithm == std::string("parmetis") ||
796  algorithm == std::string("scotch") ||
797  algorithm == std::string("ptscotch"))
798  {
799 
800  //modelType_ = GraphModelType;
801  modelAvail_[GraphModelType]=true;
802 
803  algName_ = algorithm;
804  removeSelfEdges = true;
805  needConsecutiveGlobalIds = true;
806  }
807  else if (algorithm == std::string("patoh") ||
808  algorithm == std::string("phg"))
809  {
810  // if ((modelType_ != GraphModelType) &&
811  // (modelType_ != HypergraphModelType) )
812  if ((modelAvail_[GraphModelType]==false) &&
813  (modelAvail_[HypergraphModelType]==false) )
814  {
815  //modelType_ = HypergraphModelType;
816  modelAvail_[HypergraphModelType]=true;
817  }
818  algName_ = algorithm;
819  }
820 #ifdef INCLUDE_ZOLTAN2_EXPERIMENTAL_WOLF
821  else if (algorithm == std::string("nd"))
822  {
823  modelAvail_[GraphModelType]=true;
824  modelAvail_[CoordinateModelType]=true;
825  algName_ = algorithm;
826  }
827 #endif
828  else
829  {
830  // Parameter list should ensure this does not happen.
831  throw std::logic_error("parameter list algorithm is invalid");
832  }
833  }
834  else if (model != defString)
835  {
836  // Figure out the algorithm suggested by the model.
837  if (model == std::string("hypergraph"))
838  {
839  //modelType_ = HypergraphModelType;
840  modelAvail_[HypergraphModelType]=true;
841 
842  if (problemComm_->getSize() > 1)
843  algName_ = std::string("phg");
844  else
845  algName_ = std::string("patoh");
846  }
847  else if (model == std::string("graph"))
848  {
849  //modelType_ = GraphModelType;
850  modelAvail_[GraphModelType]=true;
851 
852 #ifdef HAVE_ZOLTAN2_SCOTCH
853  if (problemComm_->getSize() > 1)
854  algName_ = std::string("ptscotch");
855  else
856  algName_ = std::string("scotch");
857  removeSelfEdges = true;
858  needConsecutiveGlobalIds = true;
859 #else
860 #ifdef HAVE_ZOLTAN2_PARMETIS
861  if (problemComm_->getSize() > 1)
862  algName_ = std::string("parmetis");
863  else
864  algName_ = std::string("metis");
865  removeSelfEdges = true;
866  needConsecutiveGlobalIds = true;
867 #else
868  if (problemComm_->getSize() > 1)
869  algName_ = std::string("phg");
870  else
871  algName_ = std::string("patoh");
872 #endif
873 #endif
874  }
875  else if (model == std::string("geometry"))
876  {
877  //modelType_ = CoordinateModelType;
878  modelAvail_[CoordinateModelType]=true;
879 
880  algName_ = std::string("multijagged");
881  }
882  else if (model == std::string("ids"))
883  {
884  //modelType_ = IdentifierModelType;
885  modelAvail_[IdentifierModelType]=true;
886 
887  algName_ = std::string("block");
888  }
889  else
890  {
891  // Parameter list should ensure this does not happen.
892  env.localBugAssertion(__FILE__, __LINE__,
893  "parameter list model type is invalid", 1, BASIC_ASSERTION);
894  }
895  }
896  else
897  {
898  // Determine an algorithm and model suggested by the input type.
899  // TODO: this is a good time to use the time vs. quality parameter
900  // in choosing an algorithm, and setting some parameters
901 
902  if (inputType_ == MatrixAdapterType)
903  {
904  //modelType_ = HypergraphModelType;
905  modelAvail_[HypergraphModelType]=true;
906 
907  if (problemComm_->getSize() > 1)
908  algName_ = std::string("phg");
909  else
910  algName_ = std::string("patoh");
911  }
912  else if (inputType_ == GraphAdapterType ||
913  inputType_ == MeshAdapterType)
914  {
915  //modelType_ = GraphModelType;
916  modelAvail_[GraphModelType]=true;
917 
918  if (problemComm_->getSize() > 1)
919  algName_ = std::string("phg");
920  else
921  algName_ = std::string("patoh");
922  }
923  else if (inputType_ == VectorAdapterType)
924  {
925  //modelType_ = CoordinateModelType;
926  modelAvail_[CoordinateModelType]=true;
927 
928  algName_ = std::string("multijagged");
929  }
930  else if (inputType_ == IdentifierAdapterType)
931  {
932  //modelType_ = IdentifierModelType;
933  modelAvail_[IdentifierModelType]=true;
934 
935  algName_ = std::string("block");
936  }
937  else{
938  // This should never happen
939  throw std::logic_error("input type is invalid");
940  }
941  }
942 
943  // Hierarchical partitioning?
944 
945  Array<int> valueList;
946  pe = pl.getEntryPtr("topology");
947 
948  if (pe){
949  valueList = pe->getValue<Array<int> >(&valueList);
950 
951  if (!Zoltan2::noValuesAreInRangeList<int>(valueList)){
952  int *n = new int [valueList.size() + 1];
953  levelNumberParts_ = arcp(n, 0, valueList.size() + 1, true);
954  int procsPerNode = 1;
955  for (int i=0; i < valueList.size(); i++){
956  levelNumberParts_[i+1] = valueList[i];
957  procsPerNode *= valueList[i];
958  }
959  // Number of parts in the first level
960  levelNumberParts_[0] = env.numProcs_ / procsPerNode;
961 
962  if (env.numProcs_ % procsPerNode > 0)
963  levelNumberParts_[0]++;
964  }
965  }
966  else{
967  levelNumberParts_.clear();
968  }
969 
970  hierarchical_ = levelNumberParts_.size() > 0;
971 
972  // Object to be partitioned? (rows, columns, etc)
973 
974  std::string objectOfInterest(defString);
975  pe = pl.getEntryPtr("objects_to_partition");
976  if (pe)
977  objectOfInterest = pe->getValue<std::string>(&objectOfInterest);
978 
980  // Set model creation flags, if any.
981 
982  this->env_->debug(DETAILED_STATUS, " models");
983  // if (modelType_ == GraphModelType)
984  if (modelAvail_[GraphModelType]==true)
985  {
986 
987  // Any parameters in the graph sublist?
988 
989  std::string symParameter(defString);
990  pe = pl.getEntryPtr("symmetrize_graph");
991  if (pe){
992  symParameter = pe->getValue<std::string>(&symParameter);
993  if (symParameter == std::string("transpose"))
994  graphFlags_.set(SYMMETRIZE_INPUT_TRANSPOSE);
995  else if (symParameter == std::string("bipartite"))
996  graphFlags_.set(SYMMETRIZE_INPUT_BIPARTITE);
997  }
998 
999  int sgParameter = 0;
1000  pe = pl.getEntryPtr("subset_graph");
1001  if (pe)
1002  sgParameter = pe->getValue<int>(&sgParameter);
1003 
1004  if (sgParameter == 1)
1005  graphFlags_.set(BUILD_SUBSET_GRAPH);
1006 
1007  // Any special behaviors required by the algorithm?
1008 
1009  if (removeSelfEdges)
1010  graphFlags_.set(REMOVE_SELF_EDGES);
1011 
1012  if (needConsecutiveGlobalIds)
1013  graphFlags_.set(GENERATE_CONSECUTIVE_IDS);
1014 
1015  // How does user input map to vertices and edges?
1016 
1017  if (inputType_ == MatrixAdapterType){
1018  if (objectOfInterest == defString ||
1019  objectOfInterest == std::string("matrix_rows") )
1020  graphFlags_.set(VERTICES_ARE_MATRIX_ROWS);
1021  else if (objectOfInterest == std::string("matrix_columns"))
1022  graphFlags_.set(VERTICES_ARE_MATRIX_COLUMNS);
1023  else if (objectOfInterest == std::string("matrix_nonzeros"))
1024  graphFlags_.set(VERTICES_ARE_MATRIX_NONZEROS);
1025  }
1026 
1027  else if (inputType_ == MeshAdapterType){
1028  if (objectOfInterest == defString ||
1029  objectOfInterest == std::string("mesh_nodes") )
1030  graphFlags_.set(VERTICES_ARE_MESH_NODES);
1031  else if (objectOfInterest == std::string("mesh_elements"))
1032  graphFlags_.set(VERTICES_ARE_MESH_ELEMENTS);
1033  }
1034  }
1035  //MMW is it ok to remove else?
1036  // else if (modelType_ == IdentifierModelType)
1037  if (modelAvail_[IdentifierModelType]==true)
1038  {
1039 
1040  // Any special behaviors required by the algorithm?
1041 
1042  }
1043  // else if (modelType_ == CoordinateModelType)
1044  if (modelAvail_[CoordinateModelType]==true)
1045  {
1046 
1047  // Any special behaviors required by the algorithm?
1048 
1049  }
1050 
1051 
1052  if (newData ||
1053  (modelAvail_[GraphModelType]!=prevModelAvail[GraphModelType]) ||
1054  (modelAvail_[HypergraphModelType]!=prevModelAvail[HypergraphModelType])||
1055  (modelAvail_[CoordinateModelType]!=prevModelAvail[CoordinateModelType])||
1056  (modelAvail_[IdentifierModelType]!=prevModelAvail[IdentifierModelType])||
1057  // (modelType_ != previousModel) ||
1058  (graphFlags_ != previousGraphModelFlags) ||
1059  (coordFlags_ != previousCoordinateModelFlags) ||
1060  (idFlags_ != previousIdentifierModelFlags))
1061  {
1062  // Create the computational model.
1063  // Models are instantiated for base input adapter types (mesh,
1064  // matrix, graph, and so on). We pass a pointer to the input
1065  // adapter, cast as the base input type.
1066 
1067  //KDD Not sure why this shadow declaration is needed
1068  //KDD Comment out for now; revisit later if problems.
1069  //KDD const Teuchos::ParameterList pl = this->envConst_->getParameters();
1070  //bool exceptionThrow = true;
1071 
1072  if(modelAvail_[GraphModelType]==true)
1073  {
1074  this->env_->debug(DETAILED_STATUS, " building graph model");
1075  this->graphModel_ = rcp(new GraphModel<base_adapter_t>(
1076  this->baseInputAdapter_, this->envConst_, problemComm_,
1077  graphFlags_));
1078 
1079  this->baseModel_ = rcp_implicit_cast<const Model<base_adapter_t> >(
1080  this->graphModel_);
1081  }
1082 
1083  if(modelAvail_[HypergraphModelType]==true)
1084  {
1085  std::cout << "Hypergraph model not implemented yet..." << std::endl;
1086  }
1087 
1088  if(modelAvail_[CoordinateModelType]==true)
1089  {
1090  this->env_->debug(DETAILED_STATUS, " building coordinate model");
1092  this->baseInputAdapter_, this->envConst_, problemComm_,
1093  coordFlags_));
1094 
1095  this->baseModel_ = rcp_implicit_cast<const Model<base_adapter_t> >(
1096  this->coordinateModel_);
1097  }
1098 
1099  if(modelAvail_[IdentifierModelType]==true)
1100  {
1101  this->env_->debug(DETAILED_STATUS, " building identifier model");
1103  this->baseInputAdapter_, this->envConst_, problemComm_,
1104  idFlags_));
1105 
1106  this->baseModel_ = rcp_implicit_cast<const Model<base_adapter_t> >(
1107  this->identifierModel_);
1108  }
1109 
1110  this->env_->memory("After creating Model");
1111  this->env_->debug(DETAILED_STATUS, "createPartitioningProblem done");
1112  }
1113 }
1114 
1115 } // namespace Zoltan2
1116 #endif
use columns as graph vertices
RCP< GraphModel< base_adapter_t > > graphModel_
#define HELLO
ArrayRCP< const MetricValues< scalar_t > > getMetrics() const
Get the array of metrics Metrics were only computed if user requested metrics with a parameter...
const scalar_t getWeightImbalance(int idx=0) const
Returns the imbalance of the solution.
Time an algorithm (or other entity) as a whole.
RCP< const base_adapter_t > baseInputAdapter_
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.
algorithm requires consecutive ids
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
PartitioningProblem(Adapter *A, ParameterList *p, RCP< const Teuchos::Comm< int > > &comm)
Constructor where Teuchos communicator is specified.
PartitioningProblem(Adapter *A, ParameterList *p)
Constructor where communicator is the Teuchos default.
model must symmetrize input
model must symmetrize input
RCP< Algorithm< Adapter > > algorithm_
RCP< const Comm< int > > comm_
void resetParameters(ParameterList *params)
Reset the list of parameters.
Defines the PartitioningSolution class.
sub-steps, each method&#39;s entry and exit
const RCP< const Environment > & getEnvironment() const
Get the current Environment. Useful for testing.
void resetParameters(ParameterList *params)
Reset the list of parameters.
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
A class that computes and returns quality metrics.
Defines the IdentifierModel interface.
A PartitioningSolution is a solution to a partitioning problem.
use nonzeros as graph vertices
Defines the PartitioningSolutionQuality class.
BaseAdapterType
An enum to identify general types of adapters.
void printMetrics(std::ostream &os) const
Print the array of metrics.
identifier data, just a list of IDs
void localBugAssertion(const char *file, int lineNum, const char *msg, bool ok, AssertionLevel level) const
Test for valid library behavior on local process only.
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...
MachineRepresentation Class Finds the coordinate of the processor. Used to find the processor coordin...
ArrayRCP< const GraphMetricValues< scalar_t > > getGraphMetrics() const
Get the array of graphMetrics Graph metrics were only computed if user requested graph metrics with a...
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
GraphModel defines the interface required for graph models.
RCP< Environment > env_
The base class for all model classes.
IdentifierModel defines the interface for all identifier models.
Defines the GraphModel interface.
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.