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