Zoltan2
Loading...
Searching...
No Matches
Zoltan2_ImbalanceMetricsUtility.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
13#ifndef ZOLTAN2_IMBALANCEMETRICSUTILITY_HPP
14#define ZOLTAN2_IMBALANCEMETRICSUTILITY_HPP
15
18
19namespace Zoltan2{
20
64template <typename scalar_t, typename lno_t, typename part_t>
66 const RCP<const Environment> &env,
67 const RCP<const Comm<int> > &comm,
68 const ArrayView<const part_t> &part,
69 int vwgtDim,
70 const ArrayView<StridedData<lno_t, scalar_t> > &vwgts,
71 multiCriteriaNorm mcNorm,
72 part_t targetNumParts,
73 part_t &numExistingParts,
74 part_t &numNonemptyParts,
75 ArrayRCP<RCP<BaseClassMetrics<scalar_t> > > &metrics,
76 ArrayRCP<scalar_t> &globalSums)
77{
78 env->debug(DETAILED_STATUS, "Entering globalSumsByPart");
80 // Initialize return values
81
82 numExistingParts = numNonemptyParts = 0;
83
84 int numMetrics = 1; // "object count"
85 if (vwgtDim) numMetrics++; // "normed weight" or "weight 0"
86 if (vwgtDim > 1) numMetrics += vwgtDim; // "weight n"
87
88 auto next = metrics.size(); // where we will start filling
89 typedef ImbalanceMetrics<scalar_t> im_t;
90 for(int n = 0; n < numMetrics; ++n) {
91 RCP<im_t> newMetric = addNewMetric<im_t, scalar_t>(env, metrics);
92 if (vwgtDim > 1) {
93 newMetric->setNorm(multiCriteriaNorm(mcNorm));
94 }
95 }
96
98 // Figure out the global number of parts in use.
99 // Verify number of vertex weights is the same everywhere.
100
101 lno_t localNumObj = part.size();
102 part_t localNum[2], globalNum[2];
103 localNum[0] = static_cast<part_t>(vwgtDim);
104 localNum[1] = 0;
105
106 for (lno_t i=0; i < localNumObj; i++)
107 if (part[i] > localNum[1]) localNum[1] = part[i];
108
109 try{
110 reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 2,
111 localNum, globalNum);
112 }
114
115 env->globalBugAssertion(__FILE__,__LINE__,
116 "inconsistent number of vertex weights",
117 globalNum[0] == localNum[0], DEBUG_MODE_ASSERTION, comm);
118
119 part_t maxPartPlusOne = globalNum[1] + 1; // Range of possible part IDs:
120 // [0,maxPartPlusOne)
121
122 part_t globalSumSize = maxPartPlusOne * numMetrics;
123
124 scalar_t * sumBuf = new scalar_t[globalSumSize];
125 env->localMemoryAssertion(__FILE__, __LINE__, globalSumSize, sumBuf);
126 globalSums = arcp(sumBuf, 0, globalSumSize);
127
129 // Calculate the local totals by part.
130
131 scalar_t *localBuf = new scalar_t[globalSumSize];
132 env->localMemoryAssertion(__FILE__, __LINE__, globalSumSize, localBuf);
133 memset(localBuf, 0, sizeof(scalar_t) * globalSumSize);
134
135 scalar_t *obj = localBuf; // # of objects
136
137 for (lno_t i=0; i < localNumObj; i++)
138 obj[part[i]]++;
139
140 if (numMetrics > 1){
141
142 scalar_t *wgt = localBuf+maxPartPlusOne; // single normed weight or weight 0
143 try{
144 normedPartWeights<scalar_t, lno_t, part_t>(env, maxPartPlusOne,
145 part, vwgts, mcNorm, wgt);
146 }
148
149 // This code assumes the solution has the part ordered the
150 // same way as the user input. (Bug 5891 is resolved.)
151 if (vwgtDim > 1){
152 wgt += maxPartPlusOne; // individual weights
153 for (int vdim = 0; vdim < vwgtDim; vdim++){
154 for (lno_t i=0; i < localNumObj; i++)
155 wgt[part[i]] += vwgts[vdim][i];
156 wgt += maxPartPlusOne;
157 }
158 }
159 }
160
161 // Metric: local sums on process
162 metrics[next]->setName("object count");
163 metrics[next]->setMetricValue("local sum", localNumObj);
164
165 next++;
166
167 if (numMetrics > 1){
168 scalar_t *wgt = localBuf+maxPartPlusOne; // single normed weight or weight 0
169 scalar_t total = 0.0;
170
171 for (int p=0; p < maxPartPlusOne; p++){
172 total += wgt[p];
173 }
174
175 if (vwgtDim == 1)
176 metrics[next]->setName("weight 0");
177 else
178 metrics[next]->setName("normed weight");
179
180 metrics[next]->setMetricValue("local sum", total);
181
182 next++;
183
184 if (vwgtDim > 1){
185 for (int vdim = 0; vdim < vwgtDim; vdim++){
186 wgt += maxPartPlusOne;
187 total = 0.0;
188 for (int p=0; p < maxPartPlusOne; p++){
189 total += wgt[p];
190 }
191
192 std::ostringstream oss;
193 oss << "weight " << vdim;
194
195 metrics[next]->setName(oss.str());
196 metrics[next]->setMetricValue("local sum", total);
197
198 next++;
199 }
200 }
201 }
202
204 // Obtain global totals by part.
205
206 try{
207 reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, globalSumSize,
208 localBuf, sumBuf);
209 }
211
212 delete [] localBuf;
213
215 // Global sum, min, max, and average over all parts
216
217 obj = sumBuf; // # of objects
218 scalar_t min=0, max=0, sum=0;
219 next = metrics.size() - numMetrics; // MDM - this is most likely temporary
220 // to preserve the format here - we are
221 // now filling a larger array so we may
222 // not have started at 0
223
224
225 ArrayView<scalar_t> objVec(obj, maxPartPlusOne);
226 getStridedStats<scalar_t>(objVec, 1, 0, min, max, sum);
227 if (maxPartPlusOne < targetNumParts)
228 min = scalar_t(0); // Some of the target parts are empty
229
230 metrics[next]->setMetricValue("global minimum", min);
231 metrics[next]->setMetricValue("global maximum", max);
232 metrics[next]->setMetricValue("global sum", sum);
233 next++;
234
235 if (numMetrics > 1){
236 scalar_t *wgt = sumBuf + maxPartPlusOne; // single normed weight or weight 0
237
238 ArrayView<scalar_t> normedWVec(wgt, maxPartPlusOne);
239 getStridedStats<scalar_t>(normedWVec, 1, 0, min, max, sum);
240 if (maxPartPlusOne < targetNumParts)
241 min = scalar_t(0); // Some of the target parts are empty
242
243 metrics[next]->setMetricValue("global minimum", min);
244 metrics[next]->setMetricValue("global maximum", max);
245 metrics[next]->setMetricValue("global sum", sum);
246 next++;
247
248 if (vwgtDim > 1){
249 for (int vdim=0; vdim < vwgtDim; vdim++){
250 wgt += maxPartPlusOne; // individual weights
251 ArrayView<scalar_t> fromVec(wgt, maxPartPlusOne);
252 getStridedStats<scalar_t>(fromVec, 1, 0, min, max, sum);
253 if (maxPartPlusOne < targetNumParts)
254 min = scalar_t(0); // Some of the target parts are empty
255
256 metrics[next]->setMetricValue("global minimum", min);
257 metrics[next]->setMetricValue("global maximum", max);
258 metrics[next]->setMetricValue("global sum", sum);
259 next++;
260 }
261 }
262 }
263
265 // How many parts do we actually have.
266
267 numExistingParts = maxPartPlusOne;
268 obj = sumBuf; // # of objects
269
270 /*for (part_t p=nparts-1; p > 0; p--){
271 if (obj[p] > 0) break;
272 numExistingParts--;
273 }*/
274
275 numNonemptyParts = numExistingParts;
276
277 for (part_t p=0; p < numExistingParts; p++)
278 if (obj[p] == 0) numNonemptyParts--;
279
280 env->debug(DETAILED_STATUS, "Exiting globalSumsByPart");
281}
282
313template <typename Adapter>
315 const RCP<const Environment> &env,
316 const RCP<const Comm<int> > &comm,
317 multiCriteriaNorm mcNorm,
318 const Adapter *ia,
319 const PartitioningSolution<Adapter> *solution,
320 const ArrayView<const typename Adapter::part_t> &partArray,
321 const RCP<const GraphModel<typename Adapter::base_adapter_t> > &graphModel,
322 typename Adapter::part_t &numExistingParts,
323 typename Adapter::part_t &numNonemptyParts,
324 ArrayRCP<RCP<BaseClassMetrics<typename Adapter::scalar_t> > > &metrics)
325{
326 env->debug(DETAILED_STATUS, "Entering objectMetrics");
327
328 typedef typename Adapter::scalar_t scalar_t;
329 typedef typename Adapter::gno_t gno_t;
330 typedef typename Adapter::lno_t lno_t;
331 typedef typename Adapter::part_t part_t;
332 typedef typename Adapter::base_adapter_t base_adapter_t;
333 typedef StridedData<lno_t, scalar_t> sdata_t;
334
335 // Local number of objects.
336
337 size_t numLocalObjects = ia->getLocalNumIDs();
338
339 // Weights, if any, for each object.
340
341 int nWeights = ia->getNumWeightsPerID();
342 int numCriteria = (nWeights > 0 ? nWeights : 1);
343 Array<sdata_t> weights(numCriteria);
344
345 if (nWeights == 0){
346 // One set of uniform weights is implied.
347 // StridedData default constructor creates length 0 strided array.
348 weights[0] = sdata_t();
349 }
350 else{
351 // whether vertex degree is ever used as vertex weight.
352 enum BaseAdapterType adapterType = ia->adapterType();
353 bool useDegreeAsWeight = false;
354 if (adapterType == GraphAdapterType) {
355 useDegreeAsWeight = reinterpret_cast<const GraphAdapter
356 <typename Adapter::user_t, typename Adapter::userCoord_t> *>(ia)->
357 useDegreeAsWeight(0);
358 } else if (adapterType == MatrixAdapterType) {
359 useDegreeAsWeight = reinterpret_cast<const MatrixAdapter
360 <typename Adapter::user_t, typename Adapter::userCoord_t> *>(ia)->
361 useDegreeAsWeight(0);
362 } else if (adapterType == MeshAdapterType) {
363 useDegreeAsWeight =
364 reinterpret_cast<const MeshAdapter<typename Adapter::user_t> *>(ia)->
365 useDegreeAsWeight(0);
366 }
367 if (useDegreeAsWeight) {
368 ArrayView<const gno_t> Ids;
369 ArrayView<sdata_t> vwgts;
370 RCP<GraphModel<base_adapter_t> > graph;
371 if (graphModel == Teuchos::null) {
372 std::bitset<NUM_MODEL_FLAGS> modelFlags;
373 const RCP<const base_adapter_t> bia =
374 rcp(dynamic_cast<const base_adapter_t *>(ia), false);
375 graph = rcp(new GraphModel<base_adapter_t>(bia,env,comm,modelFlags));
376 graph->getVertexList(Ids, vwgts);
377 } else {
378 graphModel->getVertexList(Ids, vwgts);
379 }
380 scalar_t *wgt = new scalar_t[numLocalObjects];
381 for (int i=0; i < nWeights; i++){
382 for (size_t j=0; j < numLocalObjects; j++) {
383 wgt[j] = vwgts[i][j];
384 }
385 ArrayRCP<const scalar_t> wgtArray(wgt,0,numLocalObjects,false);
386 weights[i] = sdata_t(wgtArray, 1);
387 }
388 } else {
389 for (int i=0; i < nWeights; i++){
390 const scalar_t *wgt;
391 int stride;
392 ia->getWeightsView(wgt, stride, i);
393 ArrayRCP<const scalar_t> wgtArray(wgt,0,stride*numLocalObjects,false);
394 weights[i] = sdata_t(wgtArray, stride);
395 }
396 }
397 }
398
399 // Relative part sizes, if any, assigned to the parts.
400
401 part_t targetNumParts = comm->getSize();
402 scalar_t *psizes = NULL;
403
404 ArrayRCP<ArrayRCP<scalar_t> > partSizes(numCriteria);
405 if (solution) {
406 targetNumParts = solution->getTargetGlobalNumberOfParts();
407 for (int dim=0; dim < numCriteria; dim++){
408 if (solution->criteriaHasUniformPartSizes(dim) != true){
409 psizes = new scalar_t [targetNumParts];
410 env->localMemoryAssertion(__FILE__, __LINE__, targetNumParts, psizes);
411 for (part_t i=0; i < targetNumParts; i++){
412 psizes[i] = solution->getCriteriaPartSize(dim, i);
413 }
414 partSizes[dim] = arcp(psizes, 0, targetNumParts, true);
415 }
416 }
417 }
418
420 // Get number of parts, and the number that are non-empty.
421 // Get sums per part of objects, individual weights, and normed weight sums.
422
423 ArrayRCP<scalar_t> globalSums;
424
425 int initialMetricCount = metrics.size();
426 try{
427 globalSumsByPart<scalar_t, lno_t, part_t>(env, comm,
428 partArray, nWeights, weights.view(0, numCriteria), mcNorm,
429 targetNumParts, numExistingParts, numNonemptyParts, metrics, globalSums);
430 }
432
433 int addedMetricsCount = metrics.size() - initialMetricCount;
434
436 // Compute imbalances for the object count.
437 // (Use first index of part sizes.)
438
439 int index = initialMetricCount;
440
441 scalar_t *objCount = globalSums.getRawPtr();
442 scalar_t min, max, avg;
443 psizes=NULL;
444
445 if (partSizes[0].size() > 0)
446 psizes = partSizes[0].getRawPtr();
447
448 scalar_t gsum = metrics[index]->getMetricValue("global sum");
449 computeImbalances<scalar_t, part_t>(numExistingParts, targetNumParts, psizes,
450 gsum, objCount, min, max, avg);
451
452 metrics[index]->setMetricValue("global average", gsum / targetNumParts);
453
454 metrics[index]->setMetricValue("maximum imbalance", 1.0 + max);
455 metrics[index]->setMetricValue("average imbalance", avg);
456
458 // Compute imbalances for the normed weight sum.
459
460 scalar_t *wgts = globalSums.getRawPtr() + numExistingParts;
461
462 if (addedMetricsCount > 1){
463 ++index;
464 gsum = metrics[index]->getMetricValue("global sum");
465 computeImbalances<scalar_t, part_t>(numExistingParts, targetNumParts,
466 numCriteria, partSizes.view(0, numCriteria), gsum, wgts, min, max, avg);
467
468 metrics[index]->setMetricValue("global average", gsum / targetNumParts);
469
470 metrics[index]->setMetricValue("maximum imbalance", 1.0 + max);
471 metrics[index]->setMetricValue("average imbalance", avg);
472
473 if (addedMetricsCount > 2){
474
476 // Compute imbalances for each individual weight.
477
478 ++index;
479
480 for (int vdim=0; vdim < numCriteria; vdim++){
481 wgts += numExistingParts;
482 psizes = NULL;
483
484 if (partSizes[vdim].size() > 0)
485 psizes = partSizes[vdim].getRawPtr();
486
487 gsum = metrics[index]->getMetricValue("global sum");
488 computeImbalances<scalar_t, part_t>(numExistingParts, targetNumParts,
489 psizes, gsum, wgts, min, max, avg);
490
491 metrics[index]->setMetricValue("global average", gsum / targetNumParts);
492
493 metrics[index]->setMetricValue("maximum imbalance", 1.0 + max);
494 metrics[index]->setMetricValue("average imbalance", avg);
495 index++;
496 }
497 }
498
499 }
500 env->debug(DETAILED_STATUS, "Exiting objectMetrics");
501}
502
505template <typename scalar_t, typename part_t>
507 std::ostream &os,
508 part_t targetNumParts,
509 part_t numExistingParts,
510 part_t numNonemptyParts)
511{
512 os << "Imbalance Metrics: (" << numExistingParts << " existing parts)";
513 if (numNonemptyParts < numExistingParts) {
514 os << " (" << numNonemptyParts << " of which are non-empty)";
515 }
516 os << std::endl;
517 if (targetNumParts != numExistingParts) {
518 os << "Target number of parts is " << targetNumParts << std::endl;
519 }
521}
522
525template <typename scalar_t, typename part_t>
527 std::ostream &os,
528 part_t targetNumParts,
529 part_t numExistingParts,
530 part_t numNonemptyParts,
531 const ArrayView<RCP<BaseClassMetrics<scalar_t> > > &infoList)
532{
533 printImbalanceMetricsHeader<scalar_t, part_t>(os, targetNumParts,
534 numExistingParts,
535 numNonemptyParts);
536 for (int i=0; i < infoList.size(); i++) {
537 if (infoList[i]->getName() != METRICS_UNSET_STRING) {
538 infoList[i]->printLine(os);
539 }
540 }
541 os << std::endl;
542}
543
546template <typename scalar_t, typename part_t>
548 std::ostream &os,
549 part_t targetNumParts,
550 part_t numExistingParts,
551 part_t numNonemptyParts,
552 RCP<BaseClassMetrics<scalar_t>> metricValue)
553{
554 printImbalanceMetricsHeader<scalar_t, part_t>(os, targetNumParts,
555 numExistingParts,
556 numNonemptyParts);
557 metricValue->printLine(os);
558}
559
560} //namespace Zoltan2
561
562
563#endif
#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
map_t::global_ordinal_type gno_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