Zoltan2
Zoltan2_ImbalanceMetricsUtility.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 
49 #ifndef ZOLTAN2_IMBALANCEMETRICSUTILITY_HPP
50 #define ZOLTAN2_IMBALANCEMETRICSUTILITY_HPP
51 
54 
55 namespace Zoltan2{
56 
100 template <typename scalar_t, typename lno_t, typename part_t>
102  const RCP<const Environment> &env,
103  const RCP<const Comm<int> > &comm,
104  const ArrayView<const part_t> &part,
105  int vwgtDim,
106  const ArrayView<StridedData<lno_t, scalar_t> > &vwgts,
107  multiCriteriaNorm mcNorm,
108  part_t targetNumParts,
109  part_t &numExistingParts,
110  part_t &numNonemptyParts,
111  ArrayRCP<RCP<BaseClassMetrics<scalar_t> > > &metrics,
112  ArrayRCP<scalar_t> &globalSums)
113 {
114  env->debug(DETAILED_STATUS, "Entering globalSumsByPart");
116  // Initialize return values
117 
118  numExistingParts = numNonemptyParts = 0;
119 
120  int numMetrics = 1; // "object count"
121  if (vwgtDim) numMetrics++; // "normed weight" or "weight 0"
122  if (vwgtDim > 1) numMetrics += vwgtDim; // "weight n"
123 
124  auto next = metrics.size(); // where we will start filling
125  typedef ImbalanceMetrics<scalar_t> im_t;
126  for(int n = 0; n < numMetrics; ++n) {
127  RCP<im_t> newMetric = addNewMetric<im_t, scalar_t>(env, metrics);
128  if (vwgtDim > 1) {
129  newMetric->setNorm(multiCriteriaNorm(mcNorm));
130  }
131  }
132 
134  // Figure out the global number of parts in use.
135  // Verify number of vertex weights is the same everywhere.
136 
137  lno_t localNumObj = part.size();
138  part_t localNum[2], globalNum[2];
139  localNum[0] = static_cast<part_t>(vwgtDim);
140  localNum[1] = 0;
141 
142  for (lno_t i=0; i < localNumObj; i++)
143  if (part[i] > localNum[1]) localNum[1] = part[i];
144 
145  try{
146  reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 2,
147  localNum, globalNum);
148  }
150 
151  env->globalBugAssertion(__FILE__,__LINE__,
152  "inconsistent number of vertex weights",
153  globalNum[0] == localNum[0], DEBUG_MODE_ASSERTION, comm);
154 
155  part_t maxPartPlusOne = globalNum[1] + 1; // Range of possible part IDs:
156  // [0,maxPartPlusOne)
157 
158  part_t globalSumSize = maxPartPlusOne * numMetrics;
159 
160  scalar_t * sumBuf = new scalar_t[globalSumSize];
161  env->localMemoryAssertion(__FILE__, __LINE__, globalSumSize, sumBuf);
162  globalSums = arcp(sumBuf, 0, globalSumSize);
163 
165  // Calculate the local totals by part.
166 
167  scalar_t *localBuf = new scalar_t[globalSumSize];
168  env->localMemoryAssertion(__FILE__, __LINE__, globalSumSize, localBuf);
169  memset(localBuf, 0, sizeof(scalar_t) * globalSumSize);
170 
171  scalar_t *obj = localBuf; // # of objects
172 
173  for (lno_t i=0; i < localNumObj; i++)
174  obj[part[i]]++;
175 
176  if (numMetrics > 1){
177 
178  scalar_t *wgt = localBuf+maxPartPlusOne; // single normed weight or weight 0
179  try{
180  normedPartWeights<scalar_t, lno_t, part_t>(env, maxPartPlusOne,
181  part, vwgts, mcNorm, wgt);
182  }
184 
185  // This code assumes the solution has the part ordered the
186  // same way as the user input. (Bug 5891 is resolved.)
187  if (vwgtDim > 1){
188  wgt += maxPartPlusOne; // individual weights
189  for (int vdim = 0; vdim < vwgtDim; vdim++){
190  for (lno_t i=0; i < localNumObj; i++)
191  wgt[part[i]] += vwgts[vdim][i];
192  wgt += maxPartPlusOne;
193  }
194  }
195  }
196 
197  // Metric: local sums on process
198  metrics[next]->setName("object count");
199  metrics[next]->setMetricValue("local sum", localNumObj);
200 
201  next++;
202 
203  if (numMetrics > 1){
204  scalar_t *wgt = localBuf+maxPartPlusOne; // single normed weight or weight 0
205  scalar_t total = 0.0;
206 
207  for (int p=0; p < maxPartPlusOne; p++){
208  total += wgt[p];
209  }
210 
211  if (vwgtDim == 1)
212  metrics[next]->setName("weight 0");
213  else
214  metrics[next]->setName("normed weight");
215 
216  metrics[next]->setMetricValue("local sum", total);
217 
218  next++;
219 
220  if (vwgtDim > 1){
221  for (int vdim = 0; vdim < vwgtDim; vdim++){
222  wgt += maxPartPlusOne;
223  total = 0.0;
224  for (int p=0; p < maxPartPlusOne; p++){
225  total += wgt[p];
226  }
227 
228  std::ostringstream oss;
229  oss << "weight " << vdim;
230 
231  metrics[next]->setName(oss.str());
232  metrics[next]->setMetricValue("local sum", total);
233 
234  next++;
235  }
236  }
237  }
238 
240  // Obtain global totals by part.
241 
242  try{
243  reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, globalSumSize,
244  localBuf, sumBuf);
245  }
247 
248  delete [] localBuf;
249 
251  // Global sum, min, max, and average over all parts
252 
253  obj = sumBuf; // # of objects
254  scalar_t min=0, max=0, sum=0;
255  next = metrics.size() - numMetrics; // MDM - this is most likely temporary
256  // to preserve the format here - we are
257  // now filling a larger array so we may
258  // not have started at 0
259 
260 
261  ArrayView<scalar_t> objVec(obj, maxPartPlusOne);
262  getStridedStats<scalar_t>(objVec, 1, 0, min, max, sum);
263  if (maxPartPlusOne < targetNumParts)
264  min = scalar_t(0); // Some of the target parts are empty
265 
266  metrics[next]->setMetricValue("global minimum", min);
267  metrics[next]->setMetricValue("global maximum", max);
268  metrics[next]->setMetricValue("global sum", sum);
269  next++;
270 
271  if (numMetrics > 1){
272  scalar_t *wgt = sumBuf + maxPartPlusOne; // single normed weight or weight 0
273 
274  ArrayView<scalar_t> normedWVec(wgt, maxPartPlusOne);
275  getStridedStats<scalar_t>(normedWVec, 1, 0, min, max, sum);
276  if (maxPartPlusOne < targetNumParts)
277  min = scalar_t(0); // Some of the target parts are empty
278 
279  metrics[next]->setMetricValue("global minimum", min);
280  metrics[next]->setMetricValue("global maximum", max);
281  metrics[next]->setMetricValue("global sum", sum);
282  next++;
283 
284  if (vwgtDim > 1){
285  for (int vdim=0; vdim < vwgtDim; vdim++){
286  wgt += maxPartPlusOne; // individual weights
287  ArrayView<scalar_t> fromVec(wgt, maxPartPlusOne);
288  getStridedStats<scalar_t>(fromVec, 1, 0, min, max, sum);
289  if (maxPartPlusOne < targetNumParts)
290  min = scalar_t(0); // Some of the target parts are empty
291 
292  metrics[next]->setMetricValue("global minimum", min);
293  metrics[next]->setMetricValue("global maximum", max);
294  metrics[next]->setMetricValue("global sum", sum);
295  next++;
296  }
297  }
298  }
299 
301  // How many parts do we actually have.
302 
303  numExistingParts = maxPartPlusOne;
304  obj = sumBuf; // # of objects
305 
306  /*for (part_t p=nparts-1; p > 0; p--){
307  if (obj[p] > 0) break;
308  numExistingParts--;
309  }*/
310 
311  numNonemptyParts = numExistingParts;
312 
313  for (part_t p=0; p < numExistingParts; p++)
314  if (obj[p] == 0) numNonemptyParts--;
315 
316  env->debug(DETAILED_STATUS, "Exiting globalSumsByPart");
317 }
318 
349 template <typename Adapter>
351  const RCP<const Environment> &env,
352  const RCP<const Comm<int> > &comm,
353  multiCriteriaNorm mcNorm,
354  const Adapter *ia,
355  const PartitioningSolution<Adapter> *solution,
356  const ArrayView<const typename Adapter::part_t> &partArray,
357  const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel,
358  typename Adapter::part_t &numExistingParts,
359  typename Adapter::part_t &numNonemptyParts,
360  ArrayRCP<RCP<BaseClassMetrics<typename Adapter::scalar_t> > > &metrics)
361 {
362  env->debug(DETAILED_STATUS, "Entering objectMetrics");
363 
364  typedef typename Adapter::scalar_t scalar_t;
365  typedef typename Adapter::gno_t gno_t;
366  typedef typename Adapter::lno_t lno_t;
367  typedef typename Adapter::part_t part_t;
368  typedef typename Adapter::base_adapter_t base_adapter_t;
369  typedef StridedData<lno_t, scalar_t> sdata_t;
370 
371  // Local number of objects.
372 
373  size_t numLocalObjects = ia->getLocalNumIDs();
374 
375  // Weights, if any, for each object.
376 
377  int nWeights = ia->getNumWeightsPerID();
378  int numCriteria = (nWeights > 0 ? nWeights : 1);
379  Array<sdata_t> weights(numCriteria);
380 
381  if (nWeights == 0){
382  // One set of uniform weights is implied.
383  // StridedData default constructor creates length 0 strided array.
384  weights[0] = sdata_t();
385  }
386  else{
387  // whether vertex degree is ever used as vertex weight.
388  enum BaseAdapterType adapterType = ia->adapterType();
389  bool useDegreeAsWeight = false;
390  if (adapterType == GraphAdapterType) {
391  useDegreeAsWeight = reinterpret_cast<const GraphAdapter
392  <typename Adapter::user_t, typename Adapter::userCoord_t> *>(ia)->
393  useDegreeAsWeight(0);
394  } else if (adapterType == MatrixAdapterType) {
395  useDegreeAsWeight = reinterpret_cast<const MatrixAdapter
396  <typename Adapter::user_t, typename Adapter::userCoord_t> *>(ia)->
397  useDegreeAsWeight(0);
398  } else if (adapterType == MeshAdapterType) {
399  useDegreeAsWeight =
400  reinterpret_cast<const MeshAdapter<typename Adapter::user_t> *>(ia)->
401  useDegreeAsWeight(0);
402  }
403  if (useDegreeAsWeight) {
404  ArrayView<const gno_t> Ids;
405  ArrayView<sdata_t> vwgts;
406  RCP<GraphModel<base_adapter_t> > graph;
407  if (graphModel == Teuchos::null) {
408  std::bitset<NUM_MODEL_FLAGS> modelFlags;
409  const RCP<const base_adapter_t> bia =
410  rcp(dynamic_cast<const base_adapter_t *>(ia), false);
411  graph = rcp(new GraphModel<base_adapter_t>(bia,env,comm,modelFlags));
412  graph->getVertexList(Ids, vwgts);
413  } else {
414  graphModel->getVertexList(Ids, vwgts);
415  }
416  scalar_t *wgt = new scalar_t[numLocalObjects];
417  for (int i=0; i < nWeights; i++){
418  for (size_t j=0; j < numLocalObjects; j++) {
419  wgt[j] = vwgts[i][j];
420  }
421  ArrayRCP<const scalar_t> wgtArray(wgt,0,numLocalObjects,false);
422  weights[i] = sdata_t(wgtArray, 1);
423  }
424  } else {
425  for (int i=0; i < nWeights; i++){
426  const scalar_t *wgt;
427  int stride;
428  ia->getWeightsView(wgt, stride, i);
429  ArrayRCP<const scalar_t> wgtArray(wgt,0,stride*numLocalObjects,false);
430  weights[i] = sdata_t(wgtArray, stride);
431  }
432  }
433  }
434 
435  // Relative part sizes, if any, assigned to the parts.
436 
437  part_t targetNumParts = comm->getSize();
438  scalar_t *psizes = NULL;
439 
440  ArrayRCP<ArrayRCP<scalar_t> > partSizes(numCriteria);
441  if (solution) {
442  targetNumParts = solution->getTargetGlobalNumberOfParts();
443  for (int dim=0; dim < numCriteria; dim++){
444  if (solution->criteriaHasUniformPartSizes(dim) != true){
445  psizes = new scalar_t [targetNumParts];
446  env->localMemoryAssertion(__FILE__, __LINE__, targetNumParts, psizes);
447  for (part_t i=0; i < targetNumParts; i++){
448  psizes[i] = solution->getCriteriaPartSize(dim, i);
449  }
450  partSizes[dim] = arcp(psizes, 0, targetNumParts, true);
451  }
452  }
453  }
454 
456  // Get number of parts, and the number that are non-empty.
457  // Get sums per part of objects, individual weights, and normed weight sums.
458 
459  ArrayRCP<scalar_t> globalSums;
460 
461  int initialMetricCount = metrics.size();
462  try{
463  globalSumsByPart<scalar_t, lno_t, part_t>(env, comm,
464  partArray, nWeights, weights.view(0, numCriteria), mcNorm,
465  targetNumParts, numExistingParts, numNonemptyParts, metrics, globalSums);
466  }
468 
469  int addedMetricsCount = metrics.size() - initialMetricCount;
470 
472  // Compute imbalances for the object count.
473  // (Use first index of part sizes.)
474 
475  int index = initialMetricCount;
476 
477  scalar_t *objCount = globalSums.getRawPtr();
478  scalar_t min, max, avg;
479  psizes=NULL;
480 
481  if (partSizes[0].size() > 0)
482  psizes = partSizes[0].getRawPtr();
483 
484  scalar_t gsum = metrics[index]->getMetricValue("global sum");
485  computeImbalances<scalar_t, part_t>(numExistingParts, targetNumParts, psizes,
486  gsum, objCount, min, max, avg);
487 
488  metrics[index]->setMetricValue("global average", gsum / targetNumParts);
489 
490  metrics[index]->setMetricValue("maximum imbalance", 1.0 + max);
491  metrics[index]->setMetricValue("average imbalance", avg);
492 
494  // Compute imbalances for the normed weight sum.
495 
496  scalar_t *wgts = globalSums.getRawPtr() + numExistingParts;
497 
498  if (addedMetricsCount > 1){
499  ++index;
500  gsum = metrics[index]->getMetricValue("global sum");
501  computeImbalances<scalar_t, part_t>(numExistingParts, targetNumParts,
502  numCriteria, partSizes.view(0, numCriteria), gsum, wgts, min, max, avg);
503 
504  metrics[index]->setMetricValue("global average", gsum / targetNumParts);
505 
506  metrics[index]->setMetricValue("maximum imbalance", 1.0 + max);
507  metrics[index]->setMetricValue("average imbalance", avg);
508 
509  if (addedMetricsCount > 2){
510 
512  // Compute imbalances for each individual weight.
513 
514  ++index;
515 
516  for (int vdim=0; vdim < numCriteria; vdim++){
517  wgts += numExistingParts;
518  psizes = NULL;
519 
520  if (partSizes[vdim].size() > 0)
521  psizes = partSizes[vdim].getRawPtr();
522 
523  gsum = metrics[index]->getMetricValue("global sum");
524  computeImbalances<scalar_t, part_t>(numExistingParts, targetNumParts,
525  psizes, gsum, wgts, min, max, avg);
526 
527  metrics[index]->setMetricValue("global average", gsum / targetNumParts);
528 
529  metrics[index]->setMetricValue("maximum imbalance", 1.0 + max);
530  metrics[index]->setMetricValue("average imbalance", avg);
531  index++;
532  }
533  }
534 
535  }
536  env->debug(DETAILED_STATUS, "Exiting objectMetrics");
537 }
538 
541 template <typename scalar_t, typename part_t>
543  std::ostream &os,
544  part_t targetNumParts,
545  part_t numExistingParts,
546  part_t numNonemptyParts)
547 {
548  os << "Imbalance Metrics: (" << numExistingParts << " existing parts)";
549  if (numNonemptyParts < numExistingParts) {
550  os << " (" << numNonemptyParts << " of which are non-empty)";
551  }
552  os << std::endl;
553  if (targetNumParts != numExistingParts) {
554  os << "Target number of parts is " << targetNumParts << std::endl;
555  }
557 }
558 
561 template <typename scalar_t, typename part_t>
563  std::ostream &os,
564  part_t targetNumParts,
565  part_t numExistingParts,
566  part_t numNonemptyParts,
567  const ArrayView<RCP<BaseClassMetrics<scalar_t> > > &infoList)
568 {
569  printImbalanceMetricsHeader<scalar_t, part_t>(os, targetNumParts,
570  numExistingParts,
571  numNonemptyParts);
572  for (int i=0; i < infoList.size(); i++) {
573  if (infoList[i]->getName() != METRICS_UNSET_STRING) {
574  infoList[i]->printLine(os);
575  }
576  }
577  os << std::endl;
578 }
579 
582 template <typename scalar_t, typename part_t>
584  std::ostream &os,
585  part_t targetNumParts,
586  part_t numExistingParts,
587  part_t numNonemptyParts,
588  RCP<BaseClassMetrics<scalar_t>> metricValue)
589 {
590  printImbalanceMetricsHeader<scalar_t, part_t>(os, targetNumParts,
591  numExistingParts,
592  numNonemptyParts);
593  metricValue->printLine(os);
594 }
595 
596 } //namespace Zoltan2
597 
598 
599 #endif
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:74
#define METRICS_UNSET_STRING
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
#define Z2_THROW_OUTSIDE_ERROR(env)
Throw an error returned from outside the Zoltan2 library.
GraphAdapter defines the interface for graph-based user data.
GraphModel defines the interface required for graph models.
static void printHeader(std::ostream &os)
Print a standard header.
MatrixAdapter defines the adapter interface for matrices.
MeshAdapter defines the interface for mesh input.
A PartitioningSolution is a solution to a partitioning problem.
scalar_t getCriteriaPartSize(int idx, part_t part) const
Get the size for a given weight index and a given part.
size_t getTargetGlobalNumberOfParts() const
Returns the global number of parts desired in the solution.
bool criteriaHasUniformPartSizes(int idx) const
Determine if balancing criteria has uniform part sizes. (User can specify differing part sizes....
The StridedData class manages lists of weights or coordinates.
map_t::local_ordinal_type lno_t
Definition: mapRemotes.cpp:17
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:18
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Created by mbenlioglu on Aug 31, 2020.
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.
void globalSumsByPart(const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, const ArrayView< const part_t > &part, int vwgtDim, const ArrayView< StridedData< lno_t, scalar_t > > &vwgts, multiCriteriaNorm mcNorm, part_t targetNumParts, part_t &numExistingParts, part_t &numNonemptyParts, ArrayRCP< RCP< BaseClassMetrics< scalar_t > > > &metrics, ArrayRCP< scalar_t > &globalSums)
Given the local partitioning, compute the global sums in each part.
@ DETAILED_STATUS
sub-steps, each method's entry and exit
void printImbalanceMetrics(std::ostream &os, part_t targetNumParts, part_t numExistingParts, part_t numNonemptyParts, const ArrayView< RCP< BaseClassMetrics< scalar_t > > > &infoList)
Print out list of imbalance metrics.
@ DEBUG_MODE_ASSERTION
checks for logic errors
BaseAdapterType
An enum to identify general types of adapters.
@ GraphAdapterType
graph data
@ MatrixAdapterType
matrix data
@ MeshAdapterType
mesh data
void printImbalanceMetricsHeader(std::ostream &os, part_t targetNumParts, part_t numExistingParts, part_t numNonemptyParts)
Print out header info for imbalance metrics.
multiCriteriaNorm
Enumerator used in code for multicriteria norm choice.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
SparseMatrixAdapter_t::part_t part_t