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