Zoltan2
Loading...
Searching...
No Matches
Zoltan2_EvaluatePartition.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14#ifndef ZOLTAN2_EVALUATEPARTITION_HPP
15#define ZOLTAN2_EVALUATEPARTITION_HPP
16
23
24namespace Zoltan2{
25
33template <typename Adapter>
34class EvaluatePartition : public EvaluateBaseClass<Adapter> {
35
36private:
37
38 typedef typename Adapter::base_adapter_t base_adapter_t;
39 typedef typename Adapter::lno_t lno_t;
40 typedef typename Adapter::part_t part_t;
41 typedef typename Adapter::scalar_t scalar_t;
42
43 typedef StridedData<lno_t, scalar_t> input_t;
44
45 part_t numGlobalParts_; // desired
46 part_t targetGlobalParts_; // actual
47 part_t numNonEmpty_; // of actual
48
49 typedef BaseClassMetrics<scalar_t> base_metric_type;
50 typedef ArrayRCP<RCP<base_metric_type> > base_metric_array_type;
51 base_metric_array_type metricsBase_;
52
53protected:
54 void sharedConstructor(const Adapter *ia,
55 ParameterList *p,
56 const RCP<const Comm<int> > &problemComm,
58 const RCP<const GraphModel
59 <typename Adapter::base_adapter_t> > &graphModel);
60
61
62
75 const Adapter *ia,
76 ParameterList *p,
77 const RCP<const Comm<int> > &problemComm,
79 bool force_evaluate,
80 const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel=
81 Teuchos::null):
82 numGlobalParts_(0), targetGlobalParts_(0), numNonEmpty_(0), metricsBase_() {
83 if (force_evaluate){
84 sharedConstructor(ia, p, problemComm, soln, graphModel);
85 }
86
87 }
89 const RCP<const Environment> &_env,
90 const RCP<const Comm<int> > &_problemComm,
91 const RCP<const GraphModel<typename Adapter::base_adapter_t> > &_graph,
92 const ArrayView<const typename Adapter::part_t> &_partArray,
93 typename Adapter::part_t &_numGlobalParts,
94 ArrayRCP<RCP<BaseClassMetrics<typename Adapter::scalar_t> > > &_metricsBase,
95 ArrayRCP<typename Adapter::scalar_t> &_globalSums) {
96 globalWeightedByPart <Adapter>(_env, _problemComm, _graph,
97 _partArray, _numGlobalParts, _metricsBase, _globalSums);
98 }
99public:
101
111 const Adapter *ia,
112 ParameterList *p,
114 const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel=
115 Teuchos::null):
116 numGlobalParts_(0), targetGlobalParts_(0), numNonEmpty_(0), metricsBase_()
117 {
118 Teuchos::RCP<const Comm<int> > problemComm = Tpetra::getDefaultComm();
119 sharedConstructor(ia, p, problemComm, soln, graphModel);
120 }
121
122
133 const Adapter *ia,
134 ParameterList *p,
135 const RCP<const Comm<int> > &problemComm,
137 const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel=
138 Teuchos::null):
139 numGlobalParts_(0), targetGlobalParts_(0), numNonEmpty_(0), metricsBase_()
140 {
141 sharedConstructor(ia, p, problemComm, soln, graphModel);
142 }
143
144#ifdef HAVE_ZOLTAN2_MPI
155 const Adapter *ia,
156 ParameterList *p,
157 MPI_Comm comm,
159 const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel=
160 Teuchos::null):
161 numGlobalParts_(0), targetGlobalParts_(0), numNonEmpty_(0), metricsBase_()
162 {
163 RCP<Teuchos::OpaqueWrapper<MPI_Comm> > wrapper =
164 Teuchos::opaqueWrapper(comm);
165 RCP<const Comm<int> > problemComm =
166 rcp<const Comm<int> >(new Teuchos::MpiComm<int>(wrapper));
167 sharedConstructor(ia, p, problemComm, soln, graphModel);
168 }
169#endif
170
174 // TO DO - note that with current status it probably makes more sense to
175 // break up metricsBase_ into imbalanceMetrics_ and graphMetrics_.
176 // So instead of mixing them just keep two arrays and eliminate this function.
177 // That will clean up several places in this class.
178 ArrayView<RCP<base_metric_type>> getAllMetricsOfType(
179 std::string metricType) const {
180 // find the beginning and the end of the contiguous block
181 // the list is an ArrayRCP and must preserve any ordering
182 int beginIndex = -1;
183 int sizeOfArrayView = 0;
184 for(auto n = 0; n < metricsBase_.size(); ++n) {
185 if( metricsBase_[n]->getMetricType() == metricType ) {
186 if (beginIndex == -1) {
187 beginIndex = int(n);
188 }
189 ++sizeOfArrayView;
190 }
191 }
192 if (sizeOfArrayView == 0) {
193 return ArrayView<RCP<base_metric_type> >(); // empty array view
194 }
195 return metricsBase_.view(beginIndex, sizeOfArrayView);
196 }
197
200 scalar_t getObjectCountImbalance() const {
202 if( metrics.size() <= 0 ) {
203 throw std::logic_error("getObjectCountImbalance() was called "
204 "but no metrics data was generated for " +
205 std::string(IMBALANCE_METRICS_TYPE_NAME) + "." );
206 }
207 return metrics[0]->getMetricValue("maximum imbalance");
208 }
209
216 scalar_t getNormedImbalance() const{
218 if( metrics.size() <= 0 ) {
219 throw std::logic_error("getNormedImbalance() was called "
220 "but no metrics data was generated for " +
221 std::string(IMBALANCE_METRICS_TYPE_NAME) + "." );
222 }
223 if( metrics.size() <= 1 ) {
224 throw std::logic_error("getNormedImbalance() was called "
225 "but the normed data does not exist." );
226 }
227 return metrics[1]->getMetricValue("maximum imbalance");
228 }
229
232 scalar_t getWeightImbalance(int weightIndex) const {
233 // In this case we could have
234 // Option 1
235 // object count
236 // Option 2
237 // object count
238 // weight 0
239 // Option 3
240 // object count
241 // normed imbalance
242 // weight 0
243 // weight 1
244
245 // if we have multiple weights (meaning array size if 2 or greater) than
246 // the weights begin on index 2
247 // if we have one weight 0 (option 2) then the weights begin on index 1
249 int weight0IndexStartsAtThisArrayIndex = ( metrics.size() > 2 ) ? 2 : 1;
250 int numberOfWeights = metrics.size() - weight0IndexStartsAtThisArrayIndex;
251 int indexInArray = weight0IndexStartsAtThisArrayIndex + weightIndex;
252 if( metrics.size() <= indexInArray ) {
253 throw std::logic_error("getWeightImbalance was called with weight index "+
254 std::to_string(weightIndex) +
255 " but the maximum weight available for " +
256 std::string(IMBALANCE_METRICS_TYPE_NAME) +
257 " is weight " + std::to_string(numberOfWeights-1) +
258 "." );
259 }
260 return metrics[indexInArray]->getMetricValue("maximum imbalance");
261 }
262
265 scalar_t getMaxEdgeCut() const{
266 auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
267 if( graphMetrics.size() < 1 ) {
268 throw std::logic_error("getMaxEdgeCut() was called "
269 "but no metrics data was generated for " +
270 std::string(GRAPH_METRICS_TYPE_NAME) + "." );
271 }
272 return graphMetrics[0]->getMetricValue("global maximum");
273 }
274
277 scalar_t getMaxWeightEdgeCut(int weightIndex) const{
278 auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
279 int indexInArray = weightIndex + 1; // changed this - it used to start at 0
280 if( graphMetrics.size() <= 1 ) {
281 throw std::logic_error("getMaxWeightEdgeCut was called with "
282 "weight index " + std::to_string(weightIndex) +
283 " but no weights were available for " +
284 std::string(GRAPH_METRICS_TYPE_NAME) + "." );
285 }
286 else if( graphMetrics.size() <= indexInArray ) {
287 // the size() - 2 is because weight 0 starts at array element 1
288 // (so if the array size is 2, the maximum specified weight index is
289 // weight 0 ( 2-2 = 0 )
290 throw std::logic_error("getMaxWeightEdgeCut was called with "
291 "weight index " + std::to_string(weightIndex) +
292 " but the maximum weight available for " +
293 std::string(GRAPH_METRICS_TYPE_NAME) +
294 " is weight " +
295 std::to_string(graphMetrics.size() - 2) + "." );
296 }
297 return graphMetrics[indexInArray]->getMetricValue("global maximum");
298 }
299
302 scalar_t getTotalEdgeCut() const{
303 auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
304 if( graphMetrics.size() < 1 ) {
305 throw std::logic_error("getTotalEdgeCut() was called but no metrics "
306 "data was generated for " +
307 std::string(GRAPH_METRICS_TYPE_NAME) + "." );
308 }
309 return graphMetrics[0]->getMetricValue("global sum");
310 }
311
314 scalar_t getTotalWeightEdgeCut(int weightIndex) const{
315 auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
316 int indexInArray = weightIndex + 1; // changed this; it used to start at 0
317 if( graphMetrics.size() <= 1 ) {
318 // the size() - 2 is because weight 0 starts at array element 1 (so if
319 // the array size is 2, the maximum specified weight index is
320 // weight 0 ( 2-2 = 0 )
321 throw std::logic_error("getTotalWeightEdgeCut was called with "
322 "weight index " + std::to_string(weightIndex) +
323 " but no weights were available for " +
324 std::string(GRAPH_METRICS_TYPE_NAME) + "." );
325 }
326 else if( graphMetrics.size() <= indexInArray ) {
327 throw std::logic_error("getTotalWeightEdgeCut was called with "
328 "weight index " + std::to_string(weightIndex) +
329 " but the maximum weight available for " +
330 std::string(GRAPH_METRICS_TYPE_NAME) +
331 " is weight " +
332 std::to_string(graphMetrics.size() - 2) + "." );
333 }
334 return graphMetrics[indexInArray]->getMetricValue("global sum");
335 }
336
339 scalar_t getTotalMessages() const{
340 auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
341 if( graphMetrics.size() < 1 ) {
342 throw std::logic_error("getTotalMessages() was called but no metrics "
343 "data was generated for " +
344 std::string(GRAPH_METRICS_TYPE_NAME) + "." );
345 }
346 // TODO: Would be better to avoid hard coding the array access to [1]
347 return graphMetrics[1]->getMetricValue("global sum");
348 }
349
350
353 scalar_t getMaxMessages() const{
354 auto graphMetrics = getAllMetricsOfType(GRAPH_METRICS_TYPE_NAME);
355 if( graphMetrics.size() < 1 ) {
356 throw std::logic_error("getMaxMessages() was called but no metrics "
357 "data was generated for " +
358 std::string(GRAPH_METRICS_TYPE_NAME) + "." );
359 }
360 // TODO: Would be better to avoid hard coding the array access to [1]
361 return graphMetrics[1]->getMetricValue("global maximum");
362 }
363
366 void printMetrics(std::ostream &os) const {
367 // this could be a critical decision - do we want a blank table with
368 // headers when the list is empty - for debugging that is probably better
369 // but it's very messy to have lots of empty tables in the logs
370 ArrayView<RCP<base_metric_type>> graphMetrics =
372 if (graphMetrics.size() != 0) {
373 Zoltan2::printGraphMetrics<scalar_t, part_t>(os, targetGlobalParts_,
374 numGlobalParts_, graphMetrics);
375 }
376
377 ArrayView<RCP<base_metric_type>> imbalanceMetrics =
379 if (imbalanceMetrics.size() != 0) {
380 Zoltan2::printImbalanceMetrics<scalar_t, part_t>(os, targetGlobalParts_,
381 numGlobalParts_, numNonEmpty_, imbalanceMetrics);
382 }
383 }
384};
385
386// sharedConstructor
387template <typename Adapter>
389 const Adapter *ia,
390 ParameterList *p,
391 const RCP<const Comm<int> > &comm,
393 const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel
394)
395{
396 RCP<const Comm<int> > problemComm;
397 if (comm == Teuchos::null) {
398 problemComm = Tpetra::getDefaultComm();
399 } else {
400 problemComm = comm;
401 }
402
403 RCP<Environment> env;
404
405 try{
406 env = rcp(new Environment(*p, problemComm));
407 }
409
410 env->debug(DETAILED_STATUS, std::string("Entering EvaluatePartition"));
411 env->timerStart(MACRO_TIMERS, "Computing metrics");
412
413 // Parts to which objects are assigned.
414 size_t numLocalObjects = ia->getLocalNumIDs();
415 ArrayRCP<const part_t> parts;
416
417 if (soln) {
418 // User provided a partitioning solution; use it.
419 parts = arcp(soln->getPartListView(), 0, numLocalObjects, false);
420 env->localInputAssertion(__FILE__, __LINE__, "parts not set",
421 ((numLocalObjects == 0) || soln->getPartListView()), BASIC_ASSERTION);
422 } else {
423 // User did not provide a partitioning solution;
424 // Use input adapter's partition.
425 const part_t *tmp = NULL;
426 ia->getPartsView(tmp);
427 if (tmp != NULL)
428 parts = arcp(tmp, 0, numLocalObjects, false);
429 else {
430 // User has not provided input parts in input adapter
431 part_t *procs = new part_t[numLocalObjects];
432 for (size_t i=0;i<numLocalObjects;i++) procs[i]=problemComm->getRank();
433 parts = arcp(procs, 0, numLocalObjects, true);
434 }
435 }
436 ArrayView<const part_t> partArray = parts(0, numLocalObjects);
437
438 // When we add parameters for which weights to use, we
439 // should check those here. For now we compute metrics
440 // using all weights.
441
443 const Teuchos::ParameterEntry *pe = p->getEntryPtr("partitioning_objective");
444 if (pe){
445 std::string strChoice = pe->getValue<std::string>(&strChoice);
446 if (strChoice == std::string("multicriteria_minimize_total_weight"))
448 else if (strChoice == std::string("multicriteria_minimize_maximum_weight"))
450 }
451
452 const RCP<const base_adapter_t> bia =
453 rcp(dynamic_cast<const base_adapter_t *>(ia), false);
454
455 try{
456 imbalanceMetrics<Adapter>(env, problemComm, mcnorm, ia, soln, partArray,
457 graphModel,
458 numGlobalParts_, numNonEmpty_, metricsBase_);
459 }
461
462 if (soln)
463 targetGlobalParts_ = soln->getTargetGlobalNumberOfParts();
464 else
465 targetGlobalParts_ = problemComm->getSize();
466
467 env->timerStop(MACRO_TIMERS, "Computing metrics");
468
469 BaseAdapterType inputType = bia->adapterType();
470
471 if (inputType == GraphAdapterType ||
472 inputType == MatrixAdapterType ||
473 inputType == MeshAdapterType){
474 env->timerStart(MACRO_TIMERS, "Computing graph metrics");
475 // When we add parameters for which weights to use, we
476 // should check those here. For now we compute graph metrics
477 // using all weights.
478
479 std::bitset<NUM_MODEL_FLAGS> modelFlags;
480
481 // Create a GraphModel based on input data.
482
483 RCP<const GraphModel<base_adapter_t> > graph = graphModel;
484 if (graphModel == Teuchos::null) {
485 graph = rcp(new GraphModel<base_adapter_t>(bia, env, problemComm,
486 modelFlags));
487 }
488
489 // compute weighted cuts
490 ArrayRCP<scalar_t> globalSums;
491 try {
492 this->calculate_graph_metrics(env, problemComm, graph, partArray,
493 numGlobalParts_, metricsBase_, globalSums);
494 }
496
497 env->timerStop(MACRO_TIMERS, "Computing graph metrics");
498 }
499
500 env->debug(DETAILED_STATUS, std::string("Exiting EvaluatePartition"));
501}
502
503} // namespace Zoltan2
504
505#endif
Base class for the EvaluatePartition and EvaluateOrdering classes.
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
#define GRAPH_METRICS_TYPE_NAME
#define IMBALANCE_METRICS_TYPE_NAME
Defines the PartitioningSolution class.
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
A class that computes and returns quality metrics.
scalar_t getMaxWeightEdgeCut(int weightIndex) const
getMaxWeightEdgeCuts weighted for the specified index
scalar_t getNormedImbalance() const
Return the object normed weight imbalance. Normed imbalance is only valid if there is at least 2 elem...
void printMetrics(std::ostream &os) const
Print all metrics.
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...
scalar_t getTotalMessages() const
getTotalMessages
scalar_t getWeightImbalance(int weightIndex) const
Return the imbalance for the requested weight.
scalar_t getMaxEdgeCut() const
Return the max cut for the requested weight.
ArrayView< RCP< base_metric_type > > getAllMetricsOfType(std::string metricType) const
Return the metric list for types matching the given metric type.
scalar_t getObjectCountImbalance() const
Return the object count imbalance.
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.
scalar_t getTotalWeightEdgeCut(int weightIndex) const
getTotalWeightEdgeCut weighted for the specified index
scalar_t getTotalEdgeCut() const
getTotalEdgeCut
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.
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)
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)
scalar_t getMaxMessages() const
getMaxMessages
GraphModel defines the interface required for graph models.
A PartitioningSolution is a solution to a partitioning problem.
const part_t * getPartListView() const
Returns the part list corresponding to the global ID list.
size_t getTargetGlobalNumberOfParts() const
Returns the global number of parts desired in the solution.
The StridedData class manages lists of weights or coordinates.
Created by mbenlioglu on Aug 31, 2020.
@ MACRO_TIMERS
Time an algorithm (or other entity) as a whole.
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.
@ DETAILED_STATUS
sub-steps, each method's entry and exit
@ BASIC_ASSERTION
fast typical checks for valid arguments
BaseAdapterType
An enum to identify general types of adapters.
@ GraphAdapterType
graph data
@ MatrixAdapterType
matrix data
@ MeshAdapterType
mesh data
multiCriteriaNorm
Enumerator used in code for multicriteria norm choice.