Zoltan2
Zoltan2_EvaluatePartition.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 // @HEADER
46 //
47 // ***********************************************************************
48 //
49 // Zoltan2: A package of combinatorial algorithms for scientific computing
50 // Copyright 2012 Sandia Corporation
51 //
52 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
53 // the U.S. Government retains certain rights in this software.
54 //
55 // Redistribution and use in source and binary forms, with or without
56 // modification, are permitted provided that the following conditions are
57 // met:
58 //
59 // 1. Redistributions of source code must retain the above copyright
60 // notice, this list of conditions and the following disclaimer.
61 //
62 // 2. Redistributions in binary form must reproduce the above copyright
63 // notice, this list of conditions and the following disclaimer in the
64 // documentation and/or other materials provided with the distribution.
65 //
66 // 3. Neither the name of the Corporation nor the names of the
67 // contributors may be used to endorse or promote products derived from
68 // this software without specific prior written permission.
69 //
70 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
71 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
72 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
73 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
74 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
75 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
76 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
77 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
78 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
79 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
80 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
81 //
82 // Questions? Contact Karen Devine (kddevin@sandia.gov)
83 // Erik Boman (egboman@sandia.gov)
84 // Siva Rajamanickam (srajama@sandia.gov)
85 //
86 // ***********************************************************************
87 //
88 // @HEADER
89 
94 #ifndef ZOLTAN2_EVALUATEPARTITION_HPP
95 #define ZOLTAN2_EVALUATEPARTITION_HPP
96 
97 #include <Zoltan2_GraphMetrics.hpp>
103 
104 namespace Zoltan2{
105 
113 template <typename Adapter>
114 class EvaluatePartition : public EvaluateBaseClass<Adapter> {
115 
116 private:
117 
118  typedef typename Adapter::base_adapter_t base_adapter_t;
119  typedef typename Adapter::lno_t lno_t;
120  typedef typename Adapter::part_t part_t;
121  typedef typename Adapter::scalar_t scalar_t;
122 
124 
125  part_t numGlobalParts_; // desired
126  part_t targetGlobalParts_; // actual
127  part_t numNonEmpty_; // of actual
128 
130  typedef ArrayRCP<RCP<base_metric_type> > base_metric_array_type;
131  base_metric_array_type metricsBase_;
132 
133 protected:
134  void sharedConstructor(const Adapter *ia,
135  ParameterList *p,
136  const RCP<const Comm<int> > &problemComm,
137  const PartitioningSolution<Adapter> *soln,
138  const RCP<const GraphModel
139  <typename Adapter::base_adapter_t> > &graphModel);
140 
141 
142 
155  const Adapter *ia,
156  ParameterList *p,
157  const RCP<const Comm<int> > &problemComm,
158  const PartitioningSolution<Adapter> *soln,
159  bool force_evaluate,
160  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel=
161  Teuchos::null):
162  numGlobalParts_(0), targetGlobalParts_(0), numNonEmpty_(0), metricsBase_() {
163  if (force_evaluate){
164  sharedConstructor(ia, p, problemComm, soln, graphModel);
165  }
166 
167  }
169  const RCP<const Environment> &_env,
170  const RCP<const Comm<int> > &_problemComm,
171  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &_graph,
172  const ArrayView<const typename Adapter::part_t> &_partArray,
173  typename Adapter::part_t &_numGlobalParts,
174  ArrayRCP<RCP<BaseClassMetrics<typename Adapter::scalar_t> > > &_metricsBase,
175  ArrayRCP<typename Adapter::scalar_t> &_globalSums){
176  globalWeightedCutsMessagesByPart <Adapter>(_env,
177  _problemComm, _graph, _partArray,
178  _numGlobalParts, _metricsBase,
179  _globalSums);
180  }
181 public:
182  virtual ~EvaluatePartition(){}
183 
193  const Adapter *ia,
194  ParameterList *p,
195  const PartitioningSolution<Adapter> *soln,
196  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel=
197  Teuchos::null):
198  numGlobalParts_(0), targetGlobalParts_(0), numNonEmpty_(0), metricsBase_()
199  {
200  Teuchos::RCP<const Comm<int> > problemComm = Tpetra::getDefaultComm();
201  sharedConstructor(ia, p, problemComm, soln, graphModel);
202  }
203 
204 
215  const Adapter *ia,
216  ParameterList *p,
217  const RCP<const Comm<int> > &problemComm,
218  const PartitioningSolution<Adapter> *soln,
219  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel=
220  Teuchos::null):
221  numGlobalParts_(0), targetGlobalParts_(0), numNonEmpty_(0), metricsBase_()
222  {
223  sharedConstructor(ia, p, problemComm, soln, graphModel);
224  }
225 
226 #ifdef HAVE_ZOLTAN2_MPI
227 
237  const Adapter *ia,
238  ParameterList *p,
239  MPI_Comm comm,
240  const PartitioningSolution<Adapter> *soln,
241  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel=
242  Teuchos::null):
243  numGlobalParts_(0), targetGlobalParts_(0), numNonEmpty_(0), metricsBase_()
244  {
245  RCP<Teuchos::OpaqueWrapper<MPI_Comm> > wrapper =
246  Teuchos::opaqueWrapper(comm);
247  RCP<const Comm<int> > problemComm =
248  rcp<const Comm<int> >(new Teuchos::MpiComm<int>(wrapper));
249  sharedConstructor(ia, p, problemComm, soln, graphModel);
250  }
251 #endif
252 
256  // TO DO - note that with current status it probably makes more sense to
257  // break up metricsBase_ into imbalanceMetrics_ and graphMetrics_.
258  // So instead of mixing them just keep two arrays and eliminate this function.
259  // That will clean up several places in this class.
260  ArrayView<RCP<base_metric_type>> getAllMetricsOfType(
261  std::string metricType) const {
262  // find the beginning and the end of the contiguous block
263  // the list is an ArrayRCP and must preserve any ordering
264  int beginIndex = -1;
265  int sizeOfArrayView = 0;
266  for(auto n = 0; n < metricsBase_.size(); ++n) {
267  if( metricsBase_[n]->getMetricType() == metricType ) {
268  if (beginIndex == -1) {
269  beginIndex = int(n);
270  }
271  ++sizeOfArrayView;
272  }
273  }
274  if (sizeOfArrayView == 0) {
275  return ArrayView<RCP<base_metric_type> >(); // empty array view
276  }
277  return metricsBase_.view(beginIndex, sizeOfArrayView);
278  }
279 
282  scalar_t getObjectCountImbalance() const {
284  if( metrics.size() <= 0 ) {
285  throw std::logic_error("getObjectCountImbalance() was called "
286  "but no metrics data was generated for " +
287  std::string(IMBALANCE_METRICS_TYPE_NAME) + "." );
288  }
289  return metrics[0]->getMetricValue("maximum imbalance");
290  }
291 
298  scalar_t getNormedImbalance() const{
300  if( metrics.size() <= 0 ) {
301  throw std::logic_error("getNormedImbalance() was called "
302  "but no metrics data was generated for " +
303  std::string(IMBALANCE_METRICS_TYPE_NAME) + "." );
304  }
305  if( metrics.size() <= 1 ) {
306  throw std::logic_error("getNormedImbalance() was called "
307  "but the normed data does not exist." );
308  }
309  return metrics[1]->getMetricValue("maximum imbalance");
310  }
311 
314  scalar_t getWeightImbalance(int weightIndex) const {
315  // In this case we could have
316  // Option 1
317  // object count
318  // Option 2
319  // object count
320  // weight 0
321  // Option 3
322  // object count
323  // normed imbalance
324  // weight 0
325  // weight 1
326 
327  // if we have multiple weights (meaning array size if 2 or greater) than
328  // the weights begin on index 2
329  // if we have one weight 0 (option 2) then the weights begin on index 1
331  int weight0IndexStartsAtThisArrayIndex = ( metrics.size() > 2 ) ? 2 : 1;
332  int numberOfWeights = metrics.size() - weight0IndexStartsAtThisArrayIndex;
333  int indexInArray = weight0IndexStartsAtThisArrayIndex + weightIndex;
334  if( metrics.size() <= indexInArray ) {
335  throw std::logic_error("getWeightImbalance was called with weight index "+
336  std::to_string(weightIndex) +
337  " but the maximum weight available for " +
338  std::string(IMBALANCE_METRICS_TYPE_NAME) +
339  " is weight " + std::to_string(numberOfWeights-1) +
340  "." );
341  }
342  return metrics[indexInArray]->getMetricValue("maximum imbalance");
343  }
344 
347  scalar_t getMaxEdgeCut() const{
348  auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
349  if( graphMetrics.size() < 1 ) {
350  throw std::logic_error("getMaxEdgeCut() was called "
351  "but no metrics data was generated for " +
352  std::string(GRAPH_METRICS_TYPE_NAME) + "." );
353  }
354  return graphMetrics[0]->getMetricValue("global maximum");
355  }
356 
359  scalar_t getMaxWeightEdgeCut(int weightIndex) const{
360  auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
361  int indexInArray = weightIndex + 1; // changed this - it used to start at 0
362  if( graphMetrics.size() <= 1 ) {
363  throw std::logic_error("getMaxWeightEdgeCut was called with "
364  "weight index " + std::to_string(weightIndex) +
365  " but no weights were available for " +
366  std::string(GRAPH_METRICS_TYPE_NAME) + "." );
367  }
368  else if( graphMetrics.size() <= indexInArray ) {
369  // the size() - 2 is because weight 0 starts at array element 1
370  // (so if the array size is 2, the maximum specified weight index is
371  // weight 0 ( 2-2 = 0 )
372  throw std::logic_error("getMaxWeightEdgeCut was called with "
373  "weight index " + std::to_string(weightIndex) +
374  " but the maximum weight available for " +
375  std::string(GRAPH_METRICS_TYPE_NAME) +
376  " is weight " +
377  std::to_string(graphMetrics.size() - 2) + "." );
378  }
379  return graphMetrics[indexInArray]->getMetricValue("global maximum");
380  }
381 
384  scalar_t getTotalEdgeCut() const{
385  auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
386  if( graphMetrics.size() < 1 ) {
387  throw std::logic_error("getTotalEdgeCut() was called but no metrics "
388  "data was generated for " +
389  std::string(GRAPH_METRICS_TYPE_NAME) + "." );
390  }
391  return graphMetrics[0]->getMetricValue("global sum");
392  }
393 
396  scalar_t getTotalWeightEdgeCut(int weightIndex) const{
397  auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
398  int indexInArray = weightIndex + 1; // changed this; it used to start at 0
399  if( graphMetrics.size() <= 1 ) {
400  // the size() - 2 is because weight 0 starts at array element 1 (so if
401  // the array size is 2, the maximum specified weight index is
402  // weight 0 ( 2-2 = 0 )
403  throw std::logic_error("getTotalWeightEdgeCut was called with "
404  "weight index " + std::to_string(weightIndex) +
405  " but no weights were available for " +
406  std::string(GRAPH_METRICS_TYPE_NAME) + "." );
407  }
408  else if( graphMetrics.size() <= indexInArray ) {
409  throw std::logic_error("getTotalWeightEdgeCut was called with "
410  "weight index " + std::to_string(weightIndex) +
411  " but the maximum weight available for " +
412  std::string(GRAPH_METRICS_TYPE_NAME) +
413  " is weight " +
414  std::to_string(graphMetrics.size() - 2) + "." );
415  }
416  return graphMetrics[indexInArray]->getMetricValue("global sum");
417  }
418 
421  void printMetrics(std::ostream &os) const {
422  // this could be a critical decision - do we want a blank table with
423  // headers when the list is empty - for debugging that is probably better
424  // but it's very messy to have lots of empty tables in the logs
425  ArrayView<RCP<base_metric_type>> graphMetrics =
427  if (graphMetrics.size() != 0) {
428  Zoltan2::printGraphMetrics<scalar_t, part_t>(os, targetGlobalParts_,
429  numGlobalParts_, graphMetrics);
430  }
431 
432  ArrayView<RCP<base_metric_type>> imbalanceMetrics =
434  if (imbalanceMetrics.size() != 0) {
435  Zoltan2::printImbalanceMetrics<scalar_t, part_t>(os, targetGlobalParts_,
436  numGlobalParts_, numNonEmpty_, imbalanceMetrics);
437  }
438  }
439 };
440 
441 // sharedConstructor
442 template <typename Adapter>
444  const Adapter *ia,
445  ParameterList *p,
446  const RCP<const Comm<int> > &comm,
447  const PartitioningSolution<Adapter> *soln,
448  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel
449 )
450 {
451  RCP<const Comm<int> > problemComm;
452  if (comm == Teuchos::null) {
453  problemComm = Tpetra::getDefaultComm();
454  } else {
455  problemComm = comm;
456  }
457 
458  RCP<Environment> env;
459 
460  try{
461  env = rcp(new Environment(*p, problemComm));
462  }
464 
465  env->debug(DETAILED_STATUS, std::string("Entering EvaluatePartition"));
466  env->timerStart(MACRO_TIMERS, "Computing metrics");
467 
468  // Parts to which objects are assigned.
469  size_t numLocalObjects = ia->getLocalNumIDs();
470  ArrayRCP<const part_t> parts;
471 
472  if (soln) {
473  // User provided a partitioning solution; use it.
474  parts = arcp(soln->getPartListView(), 0, numLocalObjects, false);
475  env->localInputAssertion(__FILE__, __LINE__, "parts not set",
476  ((numLocalObjects == 0) || soln->getPartListView()), BASIC_ASSERTION);
477  } else {
478  // User did not provide a partitioning solution;
479  // Use input adapter's partition.
480  const part_t *tmp = NULL;
481  ia->getPartsView(tmp);
482  if (tmp != NULL)
483  parts = arcp(tmp, 0, numLocalObjects, false);
484  else {
485  // User has not provided input parts in input adapter
486  part_t *procs = new part_t[numLocalObjects];
487  for (size_t i=0;i<numLocalObjects;i++) procs[i]=problemComm->getRank();
488  parts = arcp(procs, 0, numLocalObjects, true);
489  }
490  }
491  ArrayView<const part_t> partArray = parts(0, numLocalObjects);
492 
493  // When we add parameters for which weights to use, we
494  // should check those here. For now we compute metrics
495  // using all weights.
496 
498  const Teuchos::ParameterEntry *pe = p->getEntryPtr("partitioning_objective");
499  if (pe){
500  std::string strChoice = pe->getValue<std::string>(&strChoice);
501  if (strChoice == std::string("multicriteria_minimize_total_weight"))
502  mcnorm = normMinimizeTotalWeight;
503  else if (strChoice == std::string("multicriteria_minimize_maximum_weight"))
504  mcnorm = normMinimizeMaximumWeight;
505  }
506 
507  const RCP<const base_adapter_t> bia =
508  rcp(dynamic_cast<const base_adapter_t *>(ia), false);
509 
510  try{
511  imbalanceMetrics<Adapter>(env, problemComm, mcnorm, ia, soln, partArray,
512  graphModel,
513  numGlobalParts_, numNonEmpty_, metricsBase_);
514  }
516 
517  if (soln)
518  targetGlobalParts_ = soln->getTargetGlobalNumberOfParts();
519  else
520  targetGlobalParts_ = problemComm->getSize();
521 
522  env->timerStop(MACRO_TIMERS, "Computing metrics");
523 
524  BaseAdapterType inputType = bia->adapterType();
525 
526  if (inputType == GraphAdapterType ||
527  inputType == MatrixAdapterType ||
528  inputType == MeshAdapterType){
529  env->timerStart(MACRO_TIMERS, "Computing graph metrics");
530  // When we add parameters for which weights to use, we
531  // should check those here. For now we compute graph metrics
532  // using all weights.
533 
534  std::bitset<NUM_MODEL_FLAGS> modelFlags;
535 
536  // Create a GraphModel based on input data.
537 
538  RCP<const GraphModel<base_adapter_t> > graph = graphModel;
539  if (graphModel == Teuchos::null) {
540  graph = rcp(new GraphModel<base_adapter_t>(bia, env, problemComm,
541  modelFlags));
542  }
543 
544  // compute weighted cuts
545  ArrayRCP<scalar_t> globalSums;
546  try {
547  globalWeightedCutsByPart<Adapter>(env,
548  problemComm, graph, partArray,
549  numGlobalParts_, metricsBase_,
550  globalSums);
551  this->calculate_graph_metrics(env,
552  problemComm, graph, partArray,
553  numGlobalParts_, metricsBase_,
554  globalSums);
555  }
557 
558  env->timerStop(MACRO_TIMERS, "Computing graph metrics");
559  }
560 
561  env->debug(DETAILED_STATUS, std::string("Exiting EvaluatePartition"));
562 }
563 
564 } // namespace Zoltan2
565 
566 #endif
Zoltan2::EvaluatePartition::calculate_graph_metrics
virtual void calculate_graph_metrics(const RCP< const Environment > &_env, const RCP< const Comm< int > > &_problemComm, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &_graph, const ArrayView< const typename Adapter::part_t > &_partArray, typename Adapter::part_t &_numGlobalParts, ArrayRCP< RCP< BaseClassMetrics< typename Adapter::scalar_t > > > &_metricsBase, ArrayRCP< typename Adapter::scalar_t > &_globalSums)
Definition: Zoltan2_EvaluatePartition.hpp:168
Z2_FORWARD_EXCEPTIONS
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Definition: Zoltan2_Exceptions.hpp:106
Zoltan2::MeshAdapterType
mesh data
Definition: Zoltan2_Adapter.hpp:69
Zoltan2::multiCriteriaNorm
multiCriteriaNorm
Enumerator used in code for multicriteria norm choice.
Definition: Zoltan2_Parameters.hpp:139
IMBALANCE_METRICS_TYPE_NAME
#define IMBALANCE_METRICS_TYPE_NAME
Definition: Zoltan2_ImbalanceMetrics.hpp:56
Zoltan2::PartitioningSolution
A PartitioningSolution is a solution to a partitioning problem.
Definition: Zoltan2_PartitioningSolution.hpp:55
Zoltan2_EvaluateBaseClass.hpp
Base class for the EvaluatePartition and EvaluateOrdering classes.
Zoltan2::EvaluatePartition::~EvaluatePartition
virtual ~EvaluatePartition()
Definition: Zoltan2_EvaluatePartition.hpp:182
Zoltan2::EvaluatePartition::getWeightImbalance
scalar_t getWeightImbalance(int weightIndex) const
Return the imbalance for the requested weight.
Definition: Zoltan2_EvaluatePartition.hpp:314
Zoltan2::PartitioningSolution::getPartListView
const part_t * getPartListView() const
Returns the part list corresponding to the global ID list.
Definition: Zoltan2_PartitioningSolution.hpp:372
Zoltan2::EvaluatePartition::sharedConstructor
void sharedConstructor(const Adapter *ia, ParameterList *p, const RCP< const Comm< int > > &problemComm, const PartitioningSolution< Adapter > *soln, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &graphModel)
Definition: Zoltan2_EvaluatePartition.hpp:443
Zoltan2_TestingFramework::base_adapter_t
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Definition: Zoltan2_Typedefs.hpp:168
Zoltan2_ImbalanceMetrics.hpp
Zoltan2_GraphMetrics.hpp
Zoltan2::MACRO_TIMERS
Time an algorithm (or other entity) as a whole.
Definition: Zoltan2_Parameters.hpp:120
Zoltan2::EvaluateBaseClass
Definition: Zoltan2_EvaluateBaseClass.hpp:66
Zoltan2::BASIC_ASSERTION
fast typical checks for valid arguments
Definition: Zoltan2_Parameters.hpp:83
GRAPH_METRICS_TYPE_NAME
#define GRAPH_METRICS_TYPE_NAME
Definition: Zoltan2_GraphMetrics.hpp:56
Zoltan2::EvaluatePartition::EvaluatePartition
EvaluatePartition(const Adapter *ia, ParameterList *p, const RCP< const Comm< int > > &problemComm, const PartitioningSolution< Adapter > *soln, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &graphModel=Teuchos::null)
Constructor where Teuchos communicator is specified.
Definition: Zoltan2_EvaluatePartition.hpp:214
Zoltan2::DETAILED_STATUS
sub-steps, each method's entry and exit
Definition: Zoltan2_Parameters.hpp:101
Zoltan2::normBalanceTotalMaximum
Definition: Zoltan2_Parameters.hpp:141
part_t
SparseMatrixAdapter_t::part_t part_t
Definition: partitioningTree.cpp:74
Zoltan2_GraphMetricsUtility.hpp
Zoltan2::Environment
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
Definition: Zoltan2_Environment.hpp:83
Zoltan2::EvaluatePartition::EvaluatePartition
EvaluatePartition(const Adapter *ia, ParameterList *p, const PartitioningSolution< Adapter > *soln, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &graphModel=Teuchos::null)
Constructor where communicator is Teuchos default.
Definition: Zoltan2_EvaluatePartition.hpp:192
Zoltan2_ImbalanceMetricsUtility.hpp
Zoltan2::EvaluatePartition::getAllMetricsOfType
ArrayView< RCP< base_metric_type > > getAllMetricsOfType(std::string metricType) const
Return the metric list for types matching the given metric type.
Definition: Zoltan2_EvaluatePartition.hpp:260
Zoltan2::EvaluatePartition::printMetrics
void printMetrics(std::ostream &os) const
Print all metrics.
Definition: Zoltan2_EvaluatePartition.hpp:421
Zoltan2::normMinimizeTotalWeight
Definition: Zoltan2_Parameters.hpp:140
Zoltan2::EvaluatePartition::getMaxWeightEdgeCut
scalar_t getMaxWeightEdgeCut(int weightIndex) const
getMaxWeightEdgeCuts weighted for the specified index
Definition: Zoltan2_EvaluatePartition.hpp:359
Zoltan2::EvaluatePartition::getTotalEdgeCut
scalar_t getTotalEdgeCut() const
getTotalEdgeCut
Definition: Zoltan2_EvaluatePartition.hpp:384
Zoltan2_PartitioningSolution.hpp
Defines the PartitioningSolution class.
Zoltan2::EvaluatePartition::EvaluatePartition
EvaluatePartition(const Adapter *ia, ParameterList *p, const RCP< const Comm< int > > &problemComm, const PartitioningSolution< Adapter > *soln, bool force_evaluate, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &graphModel=Teuchos::null)
Constructor where communicator is Teuchos default, and takes another parameter whether to evaluate me...
Definition: Zoltan2_EvaluatePartition.hpp:154
Zoltan2::BaseAdapterType
BaseAdapterType
An enum to identify general types of adapters.
Definition: Zoltan2_Adapter.hpp:63
Zoltan2::imbalanceMetrics
void imbalanceMetrics(const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, multiCriteriaNorm mcNorm, const Adapter *ia, const PartitioningSolution< Adapter > *solution, const ArrayView< const typename Adapter::part_t > &partArray, const RCP< const GraphModel< typename Adapter::base_adapter_t > > &graphModel, typename Adapter::part_t &numExistingParts, typename Adapter::part_t &numNonemptyParts, ArrayRCP< RCP< BaseClassMetrics< typename Adapter::scalar_t > > > &metrics)
Compute imbalance metrics for a distribution.
Definition: Zoltan2_ImbalanceMetricsUtility.hpp:349
Zoltan2::EvaluatePartition::getTotalWeightEdgeCut
scalar_t getTotalWeightEdgeCut(int weightIndex) const
getTotalWeightEdgeCut weighted for the specified index
Definition: Zoltan2_EvaluatePartition.hpp:396
Zoltan2::EvaluatePartition::getNormedImbalance
scalar_t getNormedImbalance() const
Return the object normed weight imbalance. Normed imbalance is only valid if there is at least 2 elem...
Definition: Zoltan2_EvaluatePartition.hpp:298
Zoltan2::EvaluatePartition::getMaxEdgeCut
scalar_t getMaxEdgeCut() const
Return the max cut for the requested weight.
Definition: Zoltan2_EvaluatePartition.hpp:347
Zoltan2::GraphAdapterType
graph data
Definition: Zoltan2_Adapter.hpp:68
Zoltan2::BaseClassMetrics
Definition: Zoltan2_BaseClassMetrics.hpp:62
Zoltan2::PartitioningSolution::getTargetGlobalNumberOfParts
size_t getTargetGlobalNumberOfParts() const
Returns the global number of parts desired in the solution.
Definition: Zoltan2_PartitioningSolution.hpp:165
Zoltan2
Definition: Zoltan2_AlgSerialGreedy.hpp:56
Zoltan2::MatrixAdapterType
matrix data
Definition: Zoltan2_Adapter.hpp:67
Zoltan2::EvaluatePartition::getObjectCountImbalance
scalar_t getObjectCountImbalance() const
Return the object count imbalance.
Definition: Zoltan2_EvaluatePartition.hpp:282
Zoltan2::EvaluatePartition
A class that computes and returns quality metrics.
Definition: Zoltan2_EvaluatePartition.hpp:114
Zoltan2::normMinimizeMaximumWeight
Definition: Zoltan2_Parameters.hpp:142
Zoltan2::GraphModel
GraphModel defines the interface required for graph models.
Definition: Zoltan2_GraphModel.hpp:80
Zoltan2::StridedData
The StridedData class manages lists of weights or coordinates.
Definition: Zoltan2_StridedData.hpp:76