Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Metric.cpp
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 //
46 // Test the following:
47 // EvaluatePartition class
48 // MetricValues class
49 // Metric related namespace methods
50 
51 
53 #include <Zoltan2_TestHelpers.hpp>
56 #include <stdlib.h>
57 #include <vector>
58 
59 
60 using Teuchos::ArrayRCP;
61 using Teuchos::Array;
62 using Teuchos::RCP;
63 using Teuchos::rcp;
64 using Teuchos::arcp;
65 
66 using namespace std;
67 using std::endl;
68 using std::cout;
69 
70 template<class idInput_t>
71 void doTest(RCP<const Comm<int> > comm, int numLocalObj,
72  int nWeights, int numLocalParts, bool givePartSizes);
73 
75 
76 // for testing basic input adapter
78 
79 // for testing graph adapter
80 typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t> tcrsGraph_t;
82 
83 // creates this so we can run the test suite over BasicIdentifierAdapter
84 // and XpetraCrsGraphAdapter
85 template<class idInput_t> void runTestSuite(RCP<const Comm<int> > comm) {
86  doTest<idInput_t>(comm, 10, 0, -1, false);
87  doTest<idInput_t>(comm, 10, 0, 1, false);
88  doTest<idInput_t>(comm, 10, 0, 1, true);
89  doTest<idInput_t>(comm, 10, 1, 1, false);
90  doTest<idInput_t>(comm, 10, 1, 1, true);
91  doTest<idInput_t>(comm, 10, 2, 1, false);
92  doTest<idInput_t>(comm, 10, 2, 1, true);
93  doTest<idInput_t>(comm, 10, 1, 2, true);
94  doTest<idInput_t>(comm, 10, 1, 2, false);
95  doTest<idInput_t>(comm, 10, 1, -1, false);
96  doTest<idInput_t>(comm, 10, 1, -1, true);
97  doTest<idInput_t>(comm, 10, 2, -1, false);
98 }
99 
100 int main(int narg, char *arg[])
101 {
102  Tpetra::ScopeGuard tscope(&narg, &arg);
103  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
104 
105  int rank = comm->getRank();
106 
107  // do some tests with BasicIdentifierAdapter
108  runTestSuite<basic_idInput_t>(comm);
109 
110  // now some tests with XpetraCrsGraphAdapter
111  // Note that right now these are all going to produce the same graph
112  // metrics but could be developed further
113  runTestSuite<graph_idInput_t>(comm);
114 
115  comm->barrier();
116  if (rank==0)
117  std::cout << "PASS" << std::endl;
118 }
119 
120 // to validate the results, we call evaluate_adapter_results which is
121 // templated so it can, for example, check graph metrics only for the graph
122 // adapter. Currently both basic and graph adapter setup imbalance metrics so
123 // we also do a check on that with a universal call to
124 // evaluate_imbalance_results. Currently this needs no specialization.
125 // If we add more adapters later this could be flexible to accomodate that.
126 template<class idInput_t>
127 void evaluate_imbalance_results(RCP<const Comm<int> > comm,
128  RCP<Zoltan2::EvaluatePartition<idInput_t>> metricObject, int numLocalObj,
129  int nWeights, int original_numLocalParts, bool givePartSizes) {
130  int fail = 0;
131  int rank = comm->getRank();
132 
133  zscalar_t object_count_imbalance;
134  try{
135  object_count_imbalance = metricObject->getObjectCountImbalance();
136  if(rank == 0) {
137  cout << "Object imbalance: " << object_count_imbalance << endl;
138  }
139  }
140  catch (std::exception &e){
141  fail=1;
142  }
143  TEST_FAIL_AND_EXIT(*comm, fail==0, "getObjectCountImbalance", 1);
144 
145  if (nWeights > 0){
146  try{
147  for (int i=0; i < nWeights; i++){
148  zscalar_t imb = metricObject->getWeightImbalance(i);
149  if(rank == 0){
150  cout << "Weight " << i << " imbalance: " << imb << endl;
151  }
152  }
153  }
154  catch (std::exception &e){
155  fail=10;
156  }
157  if (!fail && nWeights > 1){
158  try{
159  zscalar_t imb = metricObject->getNormedImbalance();
160  if(rank == 0){
161  cout << "Normed weight imbalance: " << imb << endl;
162  }
163  }
164  catch (std::exception &e){
165  fail=11;
166  }
167  }
168  }
169  TEST_FAIL_AND_EXIT(*comm, fail==0, "get imbalances", 1);
170 }
171 
172 template<class idInput_t>
173 void evaluate_adapter_results(RCP<const Comm<int> > comm,
174  RCP<Zoltan2::EvaluatePartition<idInput_t>> metricObject, int numLocalObj,
175  int nWeights, int original_numLocalParts, bool givePartSizes) {
176  throw std::logic_error("evaluate_result not implemented.");
177 }
178 
179 template<>
180 void evaluate_adapter_results<graph_idInput_t>(RCP<const Comm<int> > comm,
181  RCP<Zoltan2::EvaluatePartition<graph_idInput_t>> metricObject, int numLocalObj,
182  int nWeights, int original_numLocalParts, bool givePartSizes) {
183  int fail = 0;
184  int rank = comm->getRank();
185 
186  int total_edge_cut = -1;
187  try{
188  // TODO: the unweighted getTotalEdgeCut is an integer
189  // maybe the API should be changed for this and other similar cases
190  total_edge_cut = static_cast<int>(metricObject->getTotalEdgeCut());
191  if(rank == 0){
192  cout << "Total Edge Cut: " << total_edge_cut << endl;
193  }
194  }
195  catch (std::exception &e){
196  fail=1;
197  }
198  TEST_FAIL_AND_EXIT(*comm, fail==0, "getTotalEdgeCut", 1);
199 
200  int max_edge_cut = -1;
201  try{
202  max_edge_cut = static_cast<int>(metricObject->getMaxEdgeCut());
203  if(rank == 0){
204  cout << "Max Edge Cut: " << max_edge_cut << endl;
205  }
206  }
207  catch (std::exception &e){
208  fail=1;
209  }
210  TEST_FAIL_AND_EXIT(*comm, fail==0, "getMaxEdgeCut", 1);
211 
212  int total_messages = -1;
213  try{
214  total_messages = static_cast<int>(metricObject->getTotalMessages());
215  if(rank == 0){
216  cout << "Total Messages: " << total_messages << endl;
217  }
218  }
219  catch (std::exception &e){
220  fail=1;
221  }
222  TEST_FAIL_AND_EXIT(*comm, fail==0, "getTotalMessages", 1);
223 
224  int max_messages = -1;
225  try{
226  max_messages = static_cast<int>(metricObject->getMaxMessages());
227  if(rank == 0){
228  cout << "Max Messages: " << max_messages << endl;
229  }
230  }
231  catch (std::exception &e){
232  fail=1;
233  }
234  TEST_FAIL_AND_EXIT(*comm, fail==0, "getMaxMessages", 1);
235 
236  // Now let's check our numbers.
237  // Here we do a calculation of what getTotalEdgeCut should return based on
238  // how we set things up in create_adapter.
239  // Currently the algorithm simply has every object create two links, one
240  // to the first global id and one to the last.
241  // Two of the procs will contain one of those global ids so they only have
242  // edge cuts equal to numLocalObjs to send to the other.
243  // So that is the (2 * numLocalObjs) term.
244  // All other procs will contain neither of those global ids so they have
245  // to send their objects to two procs.
246  // So that is the ((num_procs-2) * numLocalObjs * 2 term.
247  int num_procs = comm->getSize();
248  int expected_total_edge_cuts = (num_procs == 1) ? 0 :
249  (2 * numLocalObj) + ((num_procs-2) * numLocalObj * 2);
250  TEST_FAIL_AND_EXIT(*comm, total_edge_cut == expected_total_edge_cuts,
251  "getTotalEdgeCut is not the expected ", 1);
252 
253  // we can also calculate max edge cuts
254  // if num_procs 1, then it's 0
255  // if num_procs 2, then it's numLocalObjs
256  // otherwise it's 2 * numLocalObjs because at least one proc is sending
257  // to two other procs
258  int expected_max_edge_cuts = (num_procs == 1) ? 0 :
259  (num_procs == 2) ? numLocalObj : numLocalObj * 2;
260  TEST_FAIL_AND_EXIT(*comm, max_edge_cut == expected_max_edge_cuts,
261  "getMaxEdgeCut is not the expected value", 1);
262 
263  // now check total messages - in present form we can simply divide but in
264  // future things could be generalized
265  int expected_total_messages = expected_total_edge_cuts / numLocalObj;
266  TEST_FAIL_AND_EXIT(*comm, total_messages == expected_total_messages,
267  "getTotalMessages is not the expected value", 1);
268 
269  // now check max messages - in present form we can simply divide but in
270  // future things could be more generalized
271  int expected_max_messages = expected_max_edge_cuts / numLocalObj;
272  TEST_FAIL_AND_EXIT(*comm, max_messages == expected_max_messages,
273  "getMaxMessages is not the expected value", 1);
274 
275  evaluate_imbalance_results(comm, metricObject,
276  numLocalObj, nWeights, original_numLocalParts, givePartSizes);
277 }
278 
279 // for basic_idInput_t we just call the common evaluate_imbalance_results
280 // no other specialized data to consider
281 template<>
282 void evaluate_adapter_results<basic_idInput_t>(RCP<const Comm<int> > comm,
283  RCP<Zoltan2::EvaluatePartition<basic_idInput_t>> metricObject, int numLocalObj,
284  int nWeights, int original_numLocalParts, bool givePartSizes) {
285  evaluate_imbalance_results(comm, metricObject,
286  numLocalObj, nWeights, original_numLocalParts, givePartSizes);
287 }
288 
289 template<class idInput_t>
290 idInput_t * create_adapter(RCP<const Comm<int> > comm,
291  int numLocalObj, zgno_t *myGids,
292  std::vector<const zscalar_t *> & weights,
293  std::vector<int> & strides) {
294  throw std::logic_error("create_adapter not implemented.");
295 }
296 
297 template<>
299  int numLocalObj, zgno_t *myGids,
300  std::vector<const zscalar_t *> & weights,
301  std::vector<int> & strides) {
302 
303  typedef Tpetra::Map<zlno_t, zgno_t> map_t;
304  typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t> matrix_t;
305 
306  const zgno_t gNvtx = numLocalObj * comm->getSize();
307  const Teuchos::ArrayView<const zgno_t> indexList(myGids, numLocalObj);
308  Teuchos::RCP<const map_t> map = rcp(new map_t(gNvtx, indexList, 0, comm));
309 
310  // Make some stuff in the graph
311  size_t maxRowLen = 1;
312  Teuchos::RCP<matrix_t> matrix = rcp(new matrix_t(map, maxRowLen));
313 
314  // I picked this graph as a simple test case.
315  // Something we can easily calculate the final result for as we'd like to
316  // validate this but not end up rewriting the algorithm we are testing.
317  // I have each graph element create two links to the
318  // first global index and last. That means two procs will have edge cuts
319  // equal to their numLocalObj while the rest will have 2 * numLocalObj
320  //
321  // Two of the procs will have only 1 message to send
322  // The other procs will have 2 messages to send
323  // Message max is 2
324  // Message total is going to be (2)*2 + (numProcs-2)*2
325  Teuchos::Array<zgno_t> col(2);
326  Teuchos::Array<zscalar_t> val(2); val[0] = 1.; val[1] = 1.;
327  zgno_t first_id = map->getMinAllGlobalIndex();
328  zgno_t last_id = map->getMaxAllGlobalIndex();
329  for (zlno_t i = 0; i < numLocalObj; i++) {
330  zgno_t id = map->getGlobalElement(i);
331  col[0] = first_id;
332  col[1] = last_id;
333  matrix->insertGlobalValues(id, col(), val());
334  }
335 
336  matrix->fillComplete(map, map);
337 
338  size_t nVwgts = weights.size();
339  graph_idInput_t * adapter = new graph_idInput_t(matrix->getCrsGraph(), nVwgts);
340 
341  // Set the weights
342  for (size_t j = 0; j < nVwgts; j++) {
343  adapter->setWeights(weights[j], 1, j);
344  }
345 
346  return adapter;
347 }
348 
349 template<>
351  int numLocalObj, zgno_t *myGids,
352  std::vector<const zscalar_t *> & weights,
353  std::vector<int> & strides) {
354  return new basic_idInput_t(numLocalObj, myGids, weights, strides);
355 }
356 
357 // Assumes numLocalObj is the same on every process.
358 template<class idInput_t>
359 void doTest(RCP<const Comm<int> > comm, int numLocalObj,
360  int nWeights, int numLocalParts, bool givePartSizes)
361 {
362  typedef Zoltan2::EvaluatePartition<idInput_t> quality_t;
363 
364  typedef typename idInput_t::part_t part_t;
365 
366  int rank = comm->getRank();
367 
368  int original_numLocalParts = numLocalParts; // save for log and error checking
369 
370  int nprocs = comm->getSize();
371  int fail=0;
372  srand(rank+1);
373  bool testEmptyParts = (numLocalParts < 1);
374  int numGlobalParts = 0;
375 
376  if (testEmptyParts){
377  numGlobalParts = nprocs / 2;
378  if (numGlobalParts >= 1)
379  numLocalParts = (rank < numGlobalParts ? 1 : 0);
380  else{
381  numLocalParts = 1;
382  testEmptyParts = false;
383  }
384  }
385  else{
386  numGlobalParts = nprocs * numLocalParts;
387  }
388 
389  if(rank == 0) {
390  cout << endl
391  << "Test: number of weights " << nWeights
392  << ", desired number of parts " << numGlobalParts
393  << ", calculated num local parts " << numLocalParts
394  << ", original num local parts " << original_numLocalParts
395  << (givePartSizes ? ", with differing part sizes." :
396  ", with uniform part sizes.")
397  << ", Number of procs " << nprocs
398  << ", each with " << numLocalObj << " objects, part = rank." << endl;
399  }
400 
401  // An environment. This is usually created by the problem.
402 
403  Teuchos::ParameterList pl("test list");
404  pl.set("num_local_parts", numLocalParts);
405 
406  RCP<const Zoltan2::Environment> env =
407  rcp(new Zoltan2::Environment(pl, comm));
408 
409  // A simple identifier map. Usually created by the model.
410 
411  zgno_t *myGids = new zgno_t [numLocalObj];
412  for (int i=0, x=rank*numLocalObj; i < numLocalObj; i++, x++){
413  myGids[i] = x;
414  }
415 
416  // Part sizes. Usually supplied by the user to the Problem.
417  // Then the problem supplies them to the Solution.
418 
419  int partSizeDim = (givePartSizes ? (nWeights ? nWeights : 1) : 0);
420  ArrayRCP<ArrayRCP<part_t> > ids(partSizeDim);
421  ArrayRCP<ArrayRCP<zscalar_t> > sizes(partSizeDim);
422 
423  if (givePartSizes && numLocalParts > 0){
424  part_t *myParts = new part_t [numLocalParts];
425  myParts[0] = rank * numLocalParts;
426  for (int i=1; i < numLocalParts; i++)
427  myParts[i] = myParts[i-1] + 1;
428  ArrayRCP<part_t> partNums(myParts, 0, numLocalParts, true);
429 
430  zscalar_t sizeFactor = nprocs/2 - rank;
431  if (sizeFactor < 0) sizeFactor *= -1;
432  sizeFactor += 1;
433 
434  for (int dim=0; dim < partSizeDim; dim++){
435  zscalar_t *psizes = new zscalar_t [numLocalParts];
436  for (int i=0; i < numLocalParts; i++)
437  psizes[i] = sizeFactor;
438  sizes[dim] = arcp(psizes, 0, numLocalParts, true);
439  ids[dim] = partNums;
440  }
441  }
442 
443  // An input adapter with random weights. Created by the user.
444 
445  std::vector<const zscalar_t *> weights;
446  std::vector<int> strides; // default to 1
447 
448  int len = numLocalObj*nWeights;
449  ArrayRCP<zscalar_t> wgtBuf;
450  zscalar_t *wgts = NULL;
451 
452  if (len > 0){
453  wgts = new zscalar_t [len];
454  wgtBuf = arcp(wgts, 0, len, true);
455  for (int i=0; i < len; i++)
456  wgts[i] = (zscalar_t(rand()) / zscalar_t(RAND_MAX)) + 1.0;
457  }
458 
459  for (int i=0; i < nWeights; i++, wgts+=numLocalObj)
460  weights.push_back(wgts);
461 
462  idInput_t *ia = NULL;
463 
464  try {
465  ia = create_adapter<idInput_t>(comm, numLocalObj, myGids, weights, strides);
466  }
467  catch (std::exception &e){
468  fail=1;
469  }
470 
471  TEST_FAIL_AND_EXIT(*comm, fail==0, "create adapter", 1);
472 
473  // A solution (usually created by a problem)
474 
475  RCP<Zoltan2::PartitioningSolution<idInput_t> > solution;
476 
477  try{
478  if (givePartSizes)
479  solution = rcp(new Zoltan2::PartitioningSolution<idInput_t>(
480  env, comm, nWeights,
481  ids.view(0,partSizeDim), sizes.view(0,partSizeDim)));
482  else
483  solution = rcp(new Zoltan2::PartitioningSolution<idInput_t>(
484  env, comm, nWeights));
485  }
486  catch (std::exception &e){
487  fail=1;
488  }
489 
490  TEST_FAIL_AND_EXIT(*comm, fail==0, "create solution", 1);
491 
492  // Part assignment for my objects: The algorithm usually calls this.
493 
494  part_t *partNum = new part_t [numLocalObj];
495  ArrayRCP<part_t> partAssignment(partNum, 0, numLocalObj, true);
496  for (int i=0; i < numLocalObj; i++)
497  partNum[i] = rank;
498 
499  solution->setParts(partAssignment);
500 
501  // create metric object (also usually created by a problem)
502 
503  RCP<quality_t> metricObject;
504 
505  try{
506  metricObject = rcp(new quality_t(ia, &pl, comm, solution.getRawPtr()));
507  }
508  catch (std::exception &e){
509  fail=1;
510  }
511  TEST_FAIL_AND_EXIT(*comm, fail==0, "compute metrics", 1);
512 
513  try{
514  if(rank == 0){
515  metricObject->printMetrics(cout);
516  }
517  }
518  catch (std::exception &e){
519  fail=1;
520  }
521  TEST_FAIL_AND_EXIT(*comm, fail==0, "print metrics", 1);
522 
523  // will call TEST_FAIL_AND_EXIT at each internal step
524  evaluate_adapter_results<idInput_t>(comm, metricObject,
525  numLocalObj, nWeights, original_numLocalParts, givePartSizes);
526 
527  delete ia;
528 }
idInput_t * create_adapter(RCP< const Comm< int > > comm, int numLocalObj, zgno_t *myGids, std::vector< const zscalar_t * > &weights, std::vector< int > &strides)
Definition: Metric.cpp:290
basic_idInput_t * create_adapter< basic_idInput_t >(RCP< const Comm< int > > comm, int numLocalObj, zgno_t *myGids, std::vector< const zscalar_t * > &weights, std::vector< int > &strides)
Definition: Metric.cpp:350
void evaluate_adapter_results(RCP< const Comm< int > > comm, RCP< Zoltan2::EvaluatePartition< idInput_t >> metricObject, int numLocalObj, int nWeights, int original_numLocalParts, bool givePartSizes)
Definition: Metric.cpp:173
graph_idInput_t * create_adapter< graph_idInput_t >(RCP< const Comm< int > > comm, int numLocalObj, zgno_t *myGids, std::vector< const zscalar_t * > &weights, std::vector< int > &strides)
Definition: Metric.cpp:298
void setWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to weights for the primary entity type.
int main(int narg, char *arg[])
double zscalar_t
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Provides access for Zoltan2 to Xpetra::CrsGraph data.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
int zlno_t
common code used by tests
Defines XpetraCrsGraphAdapter class.
SparseMatrixAdapter_t::part_t part_t
This class represents a collection of global Identifiers and their associated weights, if any.
A PartitioningSolution is a solution to a partitioning problem.
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tcrsGraph_t
Definition: GraphModel.cpp:83
Defines the EvaluatePartition class.
void doTest(RCP< const Comm< int > > comm, int numLocalObj, int nWeights, int numLocalParts, bool givePartSizes)
Definition: Metric.cpp:359
Zoltan2::BasicIdentifierAdapter< user_t > basic_idInput_t
Definition: Metric.cpp:77
Zoltan2::XpetraCrsGraphAdapter< tcrsGraph_t, user_t > graph_idInput_t
Definition: Metric.cpp:81
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
static const std::string fail
int zgno_t
Defines the BasicIdentifierAdapter class.
A class that computes and returns quality metrics.
void evaluate_adapter_results< graph_idInput_t >(RCP< const Comm< int > > comm, RCP< Zoltan2::EvaluatePartition< graph_idInput_t >> metricObject, int numLocalObj, int nWeights, int original_numLocalParts, bool givePartSizes)
Definition: Metric.cpp:180
void runTestSuite(RCP< const Comm< int > > comm)
Definition: Metric.cpp:85
void evaluate_imbalance_results(RCP< const Comm< int > > comm, RCP< Zoltan2::EvaluatePartition< idInput_t >> metricObject, int numLocalObj, int nWeights, int original_numLocalParts, bool givePartSizes)
Definition: Metric.cpp:127
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:74
void evaluate_adapter_results< basic_idInput_t >(RCP< const Comm< int > > comm, RCP< Zoltan2::EvaluatePartition< basic_idInput_t >> metricObject, int numLocalObj, int nWeights, int original_numLocalParts, bool givePartSizes)
Definition: Metric.cpp:282