Zoltan2
Loading...
Searching...
No Matches
rcb_C.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 <Tpetra_Map.hpp>
19#include <vector>
20#include <cstdlib>
21
26int main(int argc, char *argv[])
27{
28 Tpetra::ScopeGuard tscope(&argc, &argv);
29 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
30
31 int rank = comm->getRank();
32 int nprocs = comm->getSize();
33
34 // For convenience, we'll use the Tpetra defaults for local/global ID types
35 // Users can substitute their preferred local/global ID types
36 typedef Tpetra::Map<> Map_t;
37 typedef Map_t::local_ordinal_type localId_t;
38 typedef Map_t::global_ordinal_type globalId_t;
39
40 typedef Tpetra::Details::DefaultTypes::scalar_type scalar_t;
42
43 // TODO explain
44 typedef Zoltan2::BasicVectorAdapter<myTypes> inputAdapter_t;
46 typedef inputAdapter_t::part_t part_t;
47
49 // Create input data.
50
51 size_t localCount = 40;
52 int dim = 3;
53
54 scalar_t *coords = new scalar_t [dim * localCount];
55
56 scalar_t *x = coords;
57 scalar_t *y = x + localCount;
58 scalar_t *z = y + localCount;
59
60 // Create coordinates that range from 0 to 10.0
61
62 srand(rank);
63 scalar_t scalingFactor = 10.0 / RAND_MAX;
64
65 for (size_t i=0; i < localCount*dim; i++){
66 coords[i] = scalar_t(rand()) * scalingFactor;
67 }
68
69 // Create global ids for the coordinates.
70
71 globalId_t *globalIds = new globalId_t [localCount];
72 globalId_t offset = rank * localCount;
73
74 for (size_t i=0; i < localCount; i++)
75 globalIds[i] = offset++;
76
78 // Create parameters for an RCB problem
79
80 double tolerance = 1.1;
81
82 if (rank == 0)
83 std::cout << "Imbalance tolerance is " << tolerance << std::endl;
84
85 Teuchos::ParameterList params("test params");
86 params.set("debug_level", "basic_status");
87 params.set("debug_procs", "0");
88 params.set("error_check_level", "debug_mode_assertions");
89
90 params.set("algorithm", "rcb");
91 params.set("imbalance_tolerance", tolerance );
92 params.set("num_global_parts", nprocs);
93
96 // A simple problem with no weights.
99
100 // Create a Zoltan2 input adapter for this geometry. TODO explain
101
102 inputAdapter_t *ia1 = new inputAdapter_t(localCount,globalIds,x,y,z,1,1,1);
103
104 // Create a Zoltan2 partitioning problem
105
108
109 // Solve the problem
110
111 problem1->solve();
112
113 // create metric object where communicator is Teuchos default
114
115 quality_t *metricObject1 = new quality_t(ia1, &params, //problem1->getComm(),
116 &problem1->getSolution());
117 // Check the solution.
118
119 if (rank == 0) {
120 metricObject1->printMetrics(std::cout);
121 }
122
123 if (rank == 0){
124 scalar_t imb = metricObject1->getObjectCountImbalance();
125 if (imb <= tolerance)
126 std::cout << "pass: " << imb << std::endl;
127 else
128 std::cout << "fail: " << imb << std::endl;
129 std::cout << std::endl;
130 }
131 delete metricObject1;
132
135 // Try a problem with weights
138
139 scalar_t *weights = new scalar_t [localCount];
140 for (size_t i=0; i < localCount; i++){
141 weights[i] = 1.0 + scalar_t(rank) / scalar_t(nprocs);
142 }
143
144 // Create a Zoltan2 input adapter that includes weights.
145
146 std::vector<const scalar_t *>coordVec(2);
147 std::vector<int> coordStrides(2);
148
149 coordVec[0] = x; coordStrides[0] = 1;
150 coordVec[1] = y; coordStrides[1] = 1;
151
152 std::vector<const scalar_t *>weightVec(1);
153 std::vector<int> weightStrides(1);
154
155 weightVec[0] = weights; weightStrides[0] = 1;
156
157 inputAdapter_t *ia2=new inputAdapter_t(localCount, globalIds, coordVec,
158 coordStrides,weightVec,weightStrides);
159
160 // Create a Zoltan2 partitioning problem
161
164
165 // Solve the problem
166
167 problem2->solve();
168
169 // create metric object for MPI builds
170
171#ifdef HAVE_ZOLTAN2_MPI
172 quality_t *metricObject2 = new quality_t(ia2, &params, //problem2->getComm()
173 MPI_COMM_WORLD,
174 &problem2->getSolution());
175#else
176 quality_t *metricObject2 = new quality_t(ia2, &params, problem2->getComm(),
177 &problem2->getSolution());
178#endif
179 // Check the solution.
180
181 if (rank == 0) {
182 metricObject2->printMetrics(std::cout);
183 }
184
185 if (rank == 0){
186 scalar_t imb = metricObject2->getWeightImbalance(0);
187 if (imb <= tolerance)
188 std::cout << "pass: " << imb << std::endl;
189 else
190 std::cout << "fail: " << imb << std::endl;
191 std::cout << std::endl;
192 }
193 delete metricObject2;
194
195 if (localCount > 0){
196 delete [] weights;
197 weights = NULL;
198 }
199
202 // Try a problem with multiple weights.
205
206 // Add to the parameters the multicriteria objective.
207
208 params.set("partitioning_objective", "multicriteria_minimize_total_weight");
209
210 // Create the new weights.
211
212 weights = new scalar_t [localCount*3];
213 srand(rank);
214
215 for (size_t i=0; i < localCount*3; i+=3){
216 weights[i] = 1.0 + rank / nprocs; // weight idx 1
217 weights[i+1] = rank<nprocs/2 ? 1 : 2; // weight idx 2
218 weights[i+2] = rand()/RAND_MAX +.5; // weight idx 3
219 }
220
221 // Create a Zoltan2 input adapter with these weights.
222
223 weightVec.resize(3);
224 weightStrides.resize(3);
225
226 weightVec[0] = weights; weightStrides[0] = 3;
227 weightVec[1] = weights+1; weightStrides[1] = 3;
228 weightVec[2] = weights+2; weightStrides[2] = 3;
229
230 inputAdapter_t *ia3=new inputAdapter_t(localCount, globalIds, coordVec,
231 coordStrides,weightVec,weightStrides);
232
233 // Create a Zoltan2 partitioning problem.
234
237
238 // Solve the problem
239
240 problem3->solve();
241
242 // create metric object where Teuchos communicator is specified
243
244 quality_t *metricObject3 = new quality_t(ia3, &params, problem3->getComm(),
245 &problem3->getSolution());
246 // Check the solution.
247
248 if (rank == 0) {
249 metricObject3->printMetrics(std::cout);
250 }
251
252 if (rank == 0){
253 scalar_t imb = metricObject3->getWeightImbalance(0);
254 if (imb <= tolerance)
255 std::cout << "pass: " << imb << std::endl;
256 else
257 std::cout << "fail: " << imb << std::endl;
258 std::cout << std::endl;
259 }
260 delete metricObject3;
261
263 // Try the other multicriteria objectives.
264
265 bool dataHasChanged = false; // default is true
266
267 params.set("partitioning_objective", "multicriteria_minimize_maximum_weight");
268 problem3->resetParameters(&params);
269 problem3->solve(dataHasChanged);
270
271 // Solution changed!
272
273 metricObject3 = new quality_t(ia3, &params, problem3->getComm(),
274 &problem3->getSolution());
275 if (rank == 0){
276 metricObject3->printMetrics(std::cout);
277 scalar_t imb = metricObject3->getWeightImbalance(0);
278 if (imb <= tolerance)
279 std::cout << "pass: " << imb << std::endl;
280 else
281 std::cout << "fail: " << imb << std::endl;
282 std::cout << std::endl;
283 }
284 delete metricObject3;
285
286 params.set("partitioning_objective", "multicriteria_balance_total_maximum");
287 problem3->resetParameters(&params);
288 problem3->solve(dataHasChanged);
289
290 // Solution changed!
291
292 metricObject3 = new quality_t(ia3, &params, problem3->getComm(),
293 &problem3->getSolution());
294 if (rank == 0){
295 metricObject3->printMetrics(std::cout);
296 scalar_t imb = metricObject3->getWeightImbalance(0);
297 if (imb <= tolerance)
298 std::cout << "pass: " << imb << std::endl;
299 else
300 std::cout << "fail: " << imb << std::endl;
301 std::cout << std::endl;
302 }
303 delete metricObject3;
304
305 if (localCount > 0){
306 delete [] weights;
307 weights = NULL;
308 }
309
312 // Using part sizes, ask for some parts to be empty.
315
316 // Change the number of parts to twice the number of processes to
317 // ensure that we have more than one global part.
318
319 params.set("num_global_parts", nprocs*2);
320
321 // Using the initial problem that did not have any weights, reset
322 // parameter list, and give it some part sizes.
323
324 problem1->resetParameters(&params);
325
326 part_t partIds[2];
327 scalar_t partSizes[2];
328
329 partIds[0] = rank*2; partSizes[0] = 0;
330 partIds[1] = rank*2+1; partSizes[1] = 1;
331
332 problem1->setPartSizes(2, partIds, partSizes);
333
334 // Solve the problem. The argument "dataHasChanged" indicates
335 // that we have not changed the input data, which allows the problem
336 // so skip some work when re-solving.
337
338 dataHasChanged = false;
339
340 problem1->solve(dataHasChanged);
341
342 // Obtain the solution
343
345 problem1->getSolution();
346
347 // Check it. Part sizes should all be odd.
348
349 const part_t *partAssignments = solution4.getPartListView();
350
351 int numInEmptyParts = 0;
352 for (size_t i=0; i < localCount; i++){
353 if (partAssignments[i] % 2 == 0)
354 numInEmptyParts++;
355 }
356
357 if (rank == 0)
358 std::cout << "Request that " << nprocs << " parts be empty." <<std::endl;
359
360 // Solution changed!
361
362 metricObject1 = new quality_t(ia1, &params, //problem1->getComm(),
363 &problem1->getSolution());
364 // Check the solution.
365
366 if (rank == 0) {
367 metricObject1->printMetrics(std::cout);
368 }
369
370 if (rank == 0){
371 scalar_t imb = metricObject1->getObjectCountImbalance();
372 if (imb <= tolerance)
373 std::cout << "pass: " << imb << std::endl;
374 else
375 std::cout << "fail: " << imb << std::endl;
376 std::cout << std::endl;
377 }
378 delete metricObject1;
379
380 if (coords)
381 delete [] coords;
382
383 if (globalIds)
384 delete [] globalIds;
385
386 delete problem1;
387 delete ia1;
388 delete problem2;
389 delete ia2;
390 delete problem3;
391 delete ia3;
392
393 if (rank == 0)
394 std::cout << "PASS" << std::endl;
395}
396
Defines the BasicVectorAdapter class.
Traits for application input objects.
Defines the PartitioningProblem class.
Defines the PartitioningSolution class.
int main()
A simple class that can be the User template argument for an InputAdapter.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
A class that computes and returns quality metrics.
void printMetrics(std::ostream &os) const
Print all metrics.
scalar_t getWeightImbalance(int weightIndex) const
Return the imbalance for the requested weight.
scalar_t getObjectCountImbalance() const
Return the object count imbalance.
PartitioningProblem sets up partitioning problems for the user.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
void setPartSizes(int len, part_t *partIds, scalar_t *partSizes, bool makeCopy=true)
Set or reset relative sizes for the parts that Zoltan2 will create.
A PartitioningSolution is a solution to a partitioning problem.
const part_t * getPartListView() const
Returns the part list corresponding to the global ID list.
RCP< const Comm< int > > getComm()
Return the communicator used by the problem.
void resetParameters(ParameterList *params)
Reset the list of parameters.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
SparseMatrixAdapter_t::part_t part_t
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t