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