Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
mj_int_coordinates.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 
18 #include <Zoltan2_InputTraits.hpp>
19 #include <Tpetra_Map.hpp>
20 #include <vector>
21 #include <cstdlib>
22 
23 int main(int narg, char *arg[])
24 {
25  Tpetra::ScopeGuard scope(&narg, &arg);
26  const Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
27  int rank = comm->getRank();
28  int nprocs = comm->getSize();
29  int nFail = 0;
30 
31  typedef Tpetra::Map<> Map_t;
32  typedef Map_t::local_ordinal_type localId_t;
33  typedef Map_t::global_ordinal_type globalId_t;
34  typedef int int_scalar_t; // This is the case we are testing here.
35 
37  typedef Zoltan2::BasicVectorAdapter<myTypes> inputAdapter_t;
38 
40 
42  // Create input data.
43 
44  size_t localCount = 40;
45  int dim = 3;
46 
47  // Create coordinates that range from 0 to 9
48  int_scalar_t *coords = new int_scalar_t [dim * localCount];
49  int_scalar_t *x = coords;
50  int_scalar_t *y = x + localCount;
51  int_scalar_t *z = y + localCount;
52 
53  srand(rank);
54  for (size_t i=0; i < localCount*dim; i++)
55  coords[i] = int_scalar_t(rand() % 10);
56 
57  // Create global ids for the coordinates.
58  globalId_t *globalIds = new globalId_t [localCount];
59  globalId_t offset = rank * localCount;
60  for (size_t i=0; i < localCount; i++) globalIds[i] = offset++;
61 
63  // Create parameters for an MJ problem
64 
65  Teuchos::ParameterList params("test params");
66  params.set("debug_level", "basic_status");
67  params.set("error_check_level", "debug_mode_assertions");
68 
69  params.set("algorithm", "multijagged");
70  params.set("num_global_parts", nprocs+1);
71 
73  // Test one: No weights
74 
75  inputAdapter_t *ia1 = new inputAdapter_t(localCount,globalIds,x,y,z,1,1,1);
76 
79 
80  problem1->solve();
81 
82  quality_t *metricObject1 = new quality_t(ia1, &params, comm,
83  &problem1->getSolution());
84  if (rank == 0){
85 
86  metricObject1->printMetrics(std::cout);
87 
88  double imb = metricObject1->getObjectCountImbalance();
89  if (imb <= 1.01) // Should get perfect balance
90  std::cout << "no weights -- balance satisfied: " << imb << std::endl;
91  else {
92  std::cout << "no weights -- balance failure: " << imb << std::endl;
93  nFail++;
94  }
95  std::cout << std::endl;
96  }
97  delete metricObject1;
98  delete problem1;
99  delete ia1;
100 
102  // Test two: weighted
103  // Create a Zoltan2 input adapter that includes weights.
104 
105  int_scalar_t *weights = new int_scalar_t [localCount];
106  for (size_t i=0; i < localCount; i++) weights[i] = 1 + int_scalar_t(rank);
107 
108  std::vector<const int_scalar_t *>coordVec(3);
109  std::vector<int> coordStrides(3);
110 
111  coordVec[0] = x; coordStrides[0] = 1;
112  coordVec[1] = y; coordStrides[1] = 1;
113  coordVec[2] = z; coordStrides[2] = 1;
114 
115  std::vector<const int_scalar_t *>weightVec(1);
116  std::vector<int> weightStrides(1);
117 
118  weightVec[0] = weights; weightStrides[0] = 1;
119 
120  inputAdapter_t *ia2=new inputAdapter_t(localCount, globalIds, coordVec,
121  coordStrides,weightVec,weightStrides);
122 
125 
126  problem2->solve();
127 
128  quality_t *metricObject2 = new quality_t(ia2, &params, comm,
129  &problem2->getSolution());
130  if (rank == 0){
131 
132  metricObject2->printMetrics(std::cout);
133 
134  double imb = metricObject2->getWeightImbalance(0);
135  if (imb <= 1.01)
136  std::cout << "weighted -- balance satisfied " << imb << std::endl;
137  else {
138  std::cout << "weighted -- balance failed " << imb << std::endl;
139  nFail++;
140  }
141  std::cout << std::endl;
142  }
143  delete metricObject2;
144 
145  if (weights) delete [] weights;
146  if (coords) delete [] coords;
147  if (globalIds) delete [] globalIds;
148  delete problem2;
149  delete ia2;
150 
151  if (rank == 0) {
152  if (nFail == 0) std::cout << "PASS" << std::endl;
153  else std::cout << "FAIL: " << nFail << " tests failed" << std::endl;
154  }
155 
156  return 0;
157 }
158 
void printMetrics(std::ostream &os) const
Print all metrics.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
A simple class that can be the User template argument for an InputAdapter.
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t
int main(int narg, char **arg)
Definition: coloring1.cpp:164
Defines the PartitioningSolution class.
scalar_t getWeightImbalance(int weightIndex) const
Return the imbalance for the requested weight.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
Traits for application input objects.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
A class that computes and returns quality metrics.
Defines the BasicVectorAdapter class.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
scalar_t getObjectCountImbalance() const
Return the object count imbalance.