Zoltan2
rcb_C.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 
53 #include <Zoltan2_InputTraits.hpp>
54 #include <Tpetra_Map.hpp>
55 #include <vector>
56 #include <cstdlib>
57 
58 using namespace std;
59 using std::vector;
60 using Teuchos::RCP;
61 
66 int main(int argc, char *argv[])
67 {
68 #ifdef HAVE_ZOLTAN2_MPI
69  MPI_Init(&argc, &argv);
70  int rank, nprocs;
71  MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
72  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
73 #else
74  int rank=0, nprocs=1;
75 #endif
76 
77  // For convenience, we'll use the Tpetra defaults for local/global ID types
78  // Users can substitute their preferred local/global ID types
79  typedef Tpetra::Map<> Map_t;
80  typedef Map_t::local_ordinal_type localId_t;
81  typedef Map_t::global_ordinal_type globalId_t;
82 
83  typedef double scalar_t;
85 
86  // TODO explain
87  typedef Zoltan2::BasicVectorAdapter<myTypes> inputAdapter_t;
90 
92  // Create input data.
93 
94  size_t localCount = 40;
95  int dim = 3;
96 
97  scalar_t *coords = new scalar_t [dim * localCount];
98 
99  scalar_t *x = coords;
100  scalar_t *y = x + localCount;
101  scalar_t *z = y + localCount;
102 
103  // Create coordinates that range from 0 to 10.0
104 
105  srand(rank);
106  scalar_t scalingFactor = 10.0 / RAND_MAX;
107 
108  for (size_t i=0; i < localCount*dim; i++){
109  coords[i] = scalar_t(rand()) * scalingFactor;
110  }
111 
112  // Create global ids for the coordinates.
113 
114  globalId_t *globalIds = new globalId_t [localCount];
115  globalId_t offset = rank * localCount;
116 
117  for (size_t i=0; i < localCount; i++)
118  globalIds[i] = offset++;
119 
121  // Create parameters for an RCB problem
122 
123  double tolerance = 1.1;
124 
125  if (rank == 0)
126  std::cout << "Imbalance tolerance is " << tolerance << std::endl;
127 
128  Teuchos::ParameterList params("test params");
129  params.set("debug_level", "basic_status");
130  params.set("debug_procs", "0");
131  params.set("error_check_level", "debug_mode_assertions");
132 
133  params.set("algorithm", "rcb");
134  params.set("imbalance_tolerance", tolerance );
135  params.set("num_global_parts", nprocs);
136 
139  // A simple problem with no weights.
142 
143  // Create a Zoltan2 input adapter for this geometry. TODO explain
144 
145  inputAdapter_t *ia1 = new inputAdapter_t(localCount,globalIds,x,y,z,1,1,1);
146 
147  // Create a Zoltan2 partitioning problem
148 
151 
152  // Solve the problem
153 
154  problem1->solve();
155 
156  // create metric object where communicator is Teuchos default
157 
158  quality_t *metricObject1 = new quality_t(ia1, &params, //problem1->getComm(),
159  &problem1->getSolution());
160  // Check the solution.
161 
162  if (rank == 0) {
163  metricObject1->printMetrics(cout);
164  }
165 
166  if (rank == 0){
167  scalar_t imb = metricObject1->getObjectCountImbalance();
168  if (imb <= tolerance)
169  std::cout << "pass: " << imb << std::endl;
170  else
171  std::cout << "fail: " << imb << std::endl;
172  std::cout << std::endl;
173  }
174  delete metricObject1;
175 
178  // Try a problem with weights
181 
182  scalar_t *weights = new scalar_t [localCount];
183  for (size_t i=0; i < localCount; i++){
184  weights[i] = 1.0 + scalar_t(rank) / scalar_t(nprocs);
185  }
186 
187  // Create a Zoltan2 input adapter that includes weights.
188 
189  vector<const scalar_t *>coordVec(2);
190  vector<int> coordStrides(2);
191 
192  coordVec[0] = x; coordStrides[0] = 1;
193  coordVec[1] = y; coordStrides[1] = 1;
194 
195  vector<const scalar_t *>weightVec(1);
196  vector<int> weightStrides(1);
197 
198  weightVec[0] = weights; weightStrides[0] = 1;
199 
200  inputAdapter_t *ia2=new inputAdapter_t(localCount, globalIds, coordVec,
201  coordStrides,weightVec,weightStrides);
202 
203  // Create a Zoltan2 partitioning problem
204 
207 
208  // Solve the problem
209 
210  problem2->solve();
211 
212  // create metric object for MPI builds
213 
214 #ifdef HAVE_ZOLTAN2_MPI
215  quality_t *metricObject2 = new quality_t(ia2, &params, //problem2->getComm()
216  MPI_COMM_WORLD,
217  &problem2->getSolution());
218 #else
219  quality_t *metricObject2 = new quality_t(ia2, &params, problem2->getComm(),
220  &problem2->getSolution());
221 #endif
222  // Check the solution.
223 
224  if (rank == 0) {
225  metricObject2->printMetrics(cout);
226  }
227 
228  if (rank == 0){
229  scalar_t imb = metricObject2->getWeightImbalance(0);
230  if (imb <= tolerance)
231  std::cout << "pass: " << imb << std::endl;
232  else
233  std::cout << "fail: " << imb << std::endl;
234  std::cout << std::endl;
235  }
236  delete metricObject2;
237 
238  if (localCount > 0){
239  delete [] weights;
240  weights = NULL;
241  }
242 
245  // Try a problem with multiple weights.
248 
249  // Add to the parameters the multicriteria objective.
250 
251  params.set("partitioning_objective", "multicriteria_minimize_total_weight");
252 
253  // Create the new weights.
254 
255  weights = new scalar_t [localCount*3];
256  srand(rank);
257 
258  for (size_t i=0; i < localCount*3; i+=3){
259  weights[i] = 1.0 + rank / nprocs; // weight idx 1
260  weights[i+1] = rank<nprocs/2 ? 1 : 2; // weight idx 2
261  weights[i+2] = rand()/RAND_MAX +.5; // weight idx 3
262  }
263 
264  // Create a Zoltan2 input adapter with these weights.
265 
266  weightVec.resize(3);
267  weightStrides.resize(3);
268 
269  weightVec[0] = weights; weightStrides[0] = 3;
270  weightVec[1] = weights+1; weightStrides[1] = 3;
271  weightVec[2] = weights+2; weightStrides[2] = 3;
272 
273  inputAdapter_t *ia3=new inputAdapter_t(localCount, globalIds, coordVec,
274  coordStrides,weightVec,weightStrides);
275 
276  // Create a Zoltan2 partitioning problem.
277 
280 
281  // Solve the problem
282 
283  problem3->solve();
284 
285  // create metric object where Teuchos communicator is specified
286 
287  quality_t *metricObject3 = new quality_t(ia3, &params, problem3->getComm(),
288  &problem3->getSolution());
289  // Check the solution.
290 
291  if (rank == 0) {
292  metricObject3->printMetrics(cout);
293  }
294 
295  if (rank == 0){
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 
306  // Try the other multicriteria objectives.
307 
308  bool dataHasChanged = false; // default is true
309 
310  params.set("partitioning_objective", "multicriteria_minimize_maximum_weight");
311  problem3->resetParameters(&params);
312  problem3->solve(dataHasChanged);
313 
314  // Solution changed!
315 
316  metricObject3 = new quality_t(ia3, &params, problem3->getComm(),
317  &problem3->getSolution());
318  if (rank == 0){
319  metricObject3->printMetrics(cout);
320  scalar_t imb = metricObject3->getWeightImbalance(0);
321  if (imb <= tolerance)
322  std::cout << "pass: " << imb << std::endl;
323  else
324  std::cout << "fail: " << imb << std::endl;
325  std::cout << std::endl;
326  }
327  delete metricObject3;
328 
329  params.set("partitioning_objective", "multicriteria_balance_total_maximum");
330  problem3->resetParameters(&params);
331  problem3->solve(dataHasChanged);
332 
333  // Solution changed!
334 
335  metricObject3 = new quality_t(ia3, &params, problem3->getComm(),
336  &problem3->getSolution());
337  if (rank == 0){
338  metricObject3->printMetrics(cout);
339  scalar_t imb = metricObject3->getWeightImbalance(0);
340  if (imb <= tolerance)
341  std::cout << "pass: " << imb << std::endl;
342  else
343  std::cout << "fail: " << imb << std::endl;
344  std::cout << std::endl;
345  }
346  delete metricObject3;
347 
348  if (localCount > 0){
349  delete [] weights;
350  weights = NULL;
351  }
352 
355  // Using part sizes, ask for some parts to be empty.
358 
359  // Change the number of parts to twice the number of processes to
360  // ensure that we have more than one global part.
361 
362  params.set("num_global_parts", nprocs*2);
363 
364  // Using the initial problem that did not have any weights, reset
365  // parameter list, and give it some part sizes.
366 
367  problem1->resetParameters(&params);
368 
369  part_t partIds[2];
370  scalar_t partSizes[2];
371 
372  partIds[0] = rank*2; partSizes[0] = 0;
373  partIds[1] = rank*2+1; partSizes[1] = 1;
374 
375  problem1->setPartSizes(2, partIds, partSizes);
376 
377  // Solve the problem. The argument "dataHasChanged" indicates
378  // that we have not changed the input data, which allows the problem
379  // so skip some work when re-solving.
380 
381  dataHasChanged = false;
382 
383  problem1->solve(dataHasChanged);
384 
385  // Obtain the solution
386 
388  problem1->getSolution();
389 
390  // Check it. Part sizes should all be odd.
391 
392  const part_t *partAssignments = solution4.getPartListView();
393 
394  int numInEmptyParts = 0;
395  for (size_t i=0; i < localCount; i++){
396  if (partAssignments[i] % 2 == 0)
397  numInEmptyParts++;
398  }
399 
400  if (rank == 0)
401  std::cout << "Request that " << nprocs << " parts be empty." <<std::endl;
402 
403  // Solution changed!
404 
405  metricObject1 = new quality_t(ia1, &params, //problem1->getComm(),
406  &problem1->getSolution());
407  // Check the solution.
408 
409  if (rank == 0) {
410  metricObject1->printMetrics(cout);
411  }
412 
413  if (rank == 0){
414  scalar_t imb = metricObject1->getObjectCountImbalance();
415  if (imb <= tolerance)
416  std::cout << "pass: " << imb << std::endl;
417  else
418  std::cout << "fail: " << imb << std::endl;
419  std::cout << std::endl;
420  }
421  delete metricObject1;
422 
423  if (coords)
424  delete [] coords;
425 
426  if (globalIds)
427  delete [] globalIds;
428 
429  delete problem1;
430  delete ia1;
431  delete problem2;
432  delete ia2;
433  delete problem3;
434  delete ia3;
435 
436 #ifdef HAVE_ZOLTAN2_MPI
437  MPI_Finalize();
438 #endif
439 
440  if (rank == 0)
441  std::cout << "PASS" << std::endl;
442 }
443 
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 simple class that can be the User template argument for an InputAdapter.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
void resetParameters(ParameterList *params)
Reset the list of parameters.
Defines the PartitioningSolution class.
SparseMatrixAdapter_t::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
RCP< const Comm< int > > getComm()
Return the communicator used by the problem.
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.
const part_t * getPartListView() const
Returns the part list corresponding to the global ID list.
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.
int main(int argc, char *argv[])
Definition: rcb_C.cpp:66