Zoltan2
rcbPerformanceZ1.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_TestHelpers.hpp>
54 
55 #ifdef HAVE_ZOLTAN2_ZOLTAN
56 #include <zoltan.h>
57 #include <Teuchos_CommandLineProcessor.hpp>
58 
59 #include <vector>
60 #include <ostream>
61 #include <sstream>
62 #include <string>
63 
64 #include <Zoltan2_TestHelpers.hpp>
68 #include <GeometricGenerator.hpp>
70 using namespace std;
71 using std::vector;
72 using std::cout;
73 using std::cerr;
74 using std::endl;
75 using std::bad_alloc;
76 using Teuchos::RCP;
77 using Teuchos::rcp;
78 using Teuchos::Comm;
79 using Teuchos::ArrayView;
80 using Teuchos::ArrayRCP;
81 using Teuchos::CommandLineProcessor;
82 
83 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> tMVector_t;
84 typedef Tpetra::Map<zlno_t, zgno_t, znode_t> tMap_t;
85 
86 static ArrayRCP<ArrayRCP<zscalar_t> > weights;
87 static RCP<tMVector_t> coordinates;
88 
89 
90 const char param_comment = '#';
91 
92 string trim_right_copy(
93  const string& s,
94  const string& delimiters = " \f\n\r\t\v" )
95 {
96  return s.substr( 0, s.find_last_not_of( delimiters ) + 1 );
97 }
98 
99 string trim_left_copy(
100  const string& s,
101  const string& delimiters = " \f\n\r\t\v" )
102 {
103  return s.substr( s.find_first_not_of( delimiters ) );
104 }
105 
106 string trim_copy(
107  const string& s,
108  const string& delimiters = " \f\n\r\t\v" )
109 {
110  return trim_left_copy( trim_right_copy( s, delimiters ), delimiters );
111 }
112 
113 // Zoltan1 query functions
114 bool getArgumentValue(string &argumentid, double &argumentValue, string argumentline){
115  stringstream stream(stringstream::in | stringstream::out);
116  stream << argumentline;
117  getline(stream, argumentid, '=');
118  if (stream.eof()){
119  return false;
120  }
121  stream >> argumentValue;
122  return true;
123 }
124 
125 string convert_to_string(char *args){
126  string tmp = "";
127  for(int i = 0; args[i] != 0; i++)
128  tmp += args[i];
129  return tmp;
130 }
131 int getNumObj(void *data, int *ierr)
132 {
133  *ierr = 0;
134 
135  return coordinates->getLocalLength();
136 }
137 
138 int getDim(void *data, int *ierr)
139 {
140  *ierr = 0;
141  return 3;
142 }
143 
144 void getObjList(void *data, int numGid, int numLid,
145  zgno_t * gids, zgno_t * lids,
146  int num_wgts, float *obj_wgts, int *ierr)
147 {
148  *ierr = 0;
149  size_t localLen = coordinates->getLocalLength();
150  const zgno_t *ids = coordinates->getMap()->getNodeElementList().getRawPtr();
151  zgno_t *idsNonConst = const_cast<zgno_t *>(ids);
152 
153  if (sizeof(zgno_t) == sizeof(zgno_t)){
154  memcpy(gids, idsNonConst, sizeof(zgno_t) * localLen);
155  }
156  else{
157  for (size_t i=0; i < localLen; i++)
158  gids[i] = static_cast<zgno_t>(idsNonConst[i]);
159  }
160 
161  if (num_wgts > 0){
162  float *wgts = obj_wgts;
163  for (size_t i=0; i < localLen; i++)
164  for (int w=0; w < num_wgts; w++)
165  *wgts++ = static_cast<float>(weights[w][i]);
166  }
167 }
168 
169 void getCoords(void *data, int numGid, int numLid,
170  int numObj, zgno_t * gids, zgno_t * lids,
171  int dim, double *coords, int *ierr)
172 {
173  // I know that Zoltan asks for coordinates in gid order.
174  *ierr = 0;
175  double *val = coords;
176  const zscalar_t *x = coordinates->getData(0).getRawPtr();
177  const zscalar_t *y = coordinates->getData(1).getRawPtr();
178  const zscalar_t *z = coordinates->getData(2).getRawPtr();
179  for (zlno_t i=0; i < numObj; i++){
180  *val++ = static_cast<double>(x[i]);
181  *val++ = static_cast<double>(y[i]);
182  *val++ = static_cast<double>(z[i]);
183  }
184 }
185 
186 
187 
188 enum weightTypes{
189  upDown,
190  roundRobin,
191  increasing,
192  numWeightTypes
193 };
194 
195 ArrayRCP<zscalar_t> makeWeights(
196  const RCP<const Teuchos::Comm<int> > & comm,
197  zlno_t len, weightTypes how, zscalar_t scale, int rank)
198 {
199  zscalar_t *wgts = new zscalar_t [len];
200  if (!wgts)
201  throw bad_alloc();
202 
203  ArrayRCP<zscalar_t> wgtArray(wgts, 0, len, true);
204 
205  if (how == upDown){
206  zscalar_t val = scale + rank%2;
207  for (zlno_t i=0; i < len; i++)
208  wgts[i] = val;
209  }
210  else if (how == roundRobin){
211  for (int i=0; i < 10; i++){
212  zscalar_t val = (i + 10)*scale;
213  for (int j=i; j < len; j += 10)
214  wgts[j] = val;
215  }
216  }
217  else if (how == increasing){
218  zscalar_t val = scale + rank;
219  for (zlno_t i=0; i < len; i++)
220  wgts[i] = val;
221  }
222 
223  return wgtArray;
224 }
225 
230 const RCP<tMVector_t> getMeshCoordinates(
231  const RCP<const Teuchos::Comm<int> > & comm,
232  zgno_t numGlobalCoords)
233 {
234  int rank = comm->getRank();
235  int nprocs = comm->getSize();
236 
237  double k = log(numGlobalCoords) / 3;
238  double xdimf = exp(k) + 0.5;
239  zgno_t xdim = static_cast<int>(floor(xdimf));
240  zgno_t ydim = xdim;
241  zgno_t zdim = numGlobalCoords / (xdim*ydim);
242  zgno_t num=xdim*ydim*zdim;
243  zgno_t diff = numGlobalCoords - num;
244  zgno_t newdiff = 0;
245 
246  while (diff > 0){
247  if (zdim > xdim && zdim > ydim){
248  zdim++;
249  newdiff = diff - (xdim*ydim);
250  if (newdiff < 0)
251  if (diff < -newdiff)
252  zdim--;
253  }
254  else if (ydim > xdim && ydim > zdim){
255  ydim++;
256  newdiff = diff - (xdim*zdim);
257  if (newdiff < 0)
258  if (diff < -newdiff)
259  ydim--;
260  }
261  else{
262  xdim++;
263  newdiff = diff - (ydim*zdim);
264  if (newdiff < 0)
265  if (diff < -newdiff)
266  xdim--;
267  }
268 
269  diff = newdiff;
270  }
271 
272  num=xdim*ydim*zdim;
273  diff = numGlobalCoords - num;
274  if (diff < 0)
275  diff /= -numGlobalCoords;
276  else
277  diff /= numGlobalCoords;
278 
279  if (rank == 0){
280  if (diff > .01)
281  cout << "Warning: Difference " << diff*100 << " percent" << endl;
282  cout << "Mesh size: " << xdim << "x" << ydim << "x" <<
283  zdim << ", " << num << " vertices." << endl;
284  }
285 
286  // Divide coordinates.
287 
288  zgno_t numLocalCoords = num / nprocs;
289  zgno_t leftOver = num % nprocs;
290  zgno_t gid0 = 0;
291 
292  if (rank <= leftOver)
293  gid0 = zgno_t(rank) * (numLocalCoords+1);
294  else
295  gid0 = (leftOver * (numLocalCoords+1)) +
296  ((zgno_t(rank) - leftOver) * numLocalCoords);
297 
298  if (rank < leftOver)
299  numLocalCoords++;
300 
301  zgno_t gid1 = gid0 + numLocalCoords;
302 
303  zgno_t *ids = new zgno_t [numLocalCoords];
304  if (!ids)
305  throw bad_alloc();
306  ArrayRCP<zgno_t> idArray(ids, 0, numLocalCoords, true);
307 
308  for (zgno_t i=gid0; i < gid1; i++)
309  *ids++ = i;
310 
311  RCP<const tMap_t> idMap = rcp(
312  new tMap_t(num, idArray.view(0, numLocalCoords), 0, comm));
313 
314  // Create a Tpetra::MultiVector of coordinates.
315 
316  zscalar_t *x = new zscalar_t [numLocalCoords*3];
317  if (!x)
318  throw bad_alloc();
319  ArrayRCP<zscalar_t> coordArray(x, 0, numLocalCoords*3, true);
320 
321  zscalar_t *y = x + numLocalCoords;
322  zscalar_t *z = y + numLocalCoords;
323 
324  zgno_t xStart = 0;
325  zgno_t yStart = 0;
326  zgno_t xyPlane = xdim*ydim;
327  zgno_t zStart = gid0 / xyPlane;
328  zgno_t rem = gid0 % xyPlane;
329  if (rem > 0){
330  yStart = rem / xdim;
331  xStart = rem % xdim;
332  }
333 
334  zlno_t next = 0;
335  for (zscalar_t zval=zStart; next < numLocalCoords && zval < zdim; zval++){
336  for (zscalar_t yval=yStart; next < numLocalCoords && yval < ydim; yval++){
337  for (zscalar_t xval=xStart; next < numLocalCoords && xval < xdim; xval++){
338  x[next] = xval;
339  y[next] = yval;
340  z[next] = zval;
341  next++;
342  }
343  xStart = 0;
344  }
345  yStart = 0;
346  }
347 
348  ArrayView<const zscalar_t> xArray(x, numLocalCoords);
349  ArrayView<const zscalar_t> yArray(y, numLocalCoords);
350  ArrayView<const zscalar_t> zArray(z, numLocalCoords);
351  ArrayRCP<ArrayView<const zscalar_t> > coordinates =
352  arcp(new ArrayView<const zscalar_t> [3], 0, 3);
353  coordinates[0] = xArray;
354  coordinates[1] = yArray;
355  coordinates[2] = zArray;
356 
357  ArrayRCP<const ArrayView<const zscalar_t> > constCoords =
358  coordinates.getConst();
359 
360  RCP<tMVector_t> meshCoords = rcp(new tMVector_t(
361  idMap, constCoords.view(0,3), 3));
362 
363  return meshCoords;
364 }
365 
366 
367 void getArgVals(int argc, char **argv, int &numParts,
368  std::string &paramFile){
369 
370  for(int i = 0; i < argc; ++i){
371  string tmp = convert_to_string(argv[i]);
372  string identifier = "";
373  long long int value = -1; double fval = -1;
374  if(!getArgumentValue(identifier, fval, tmp)) continue;
375  value = (long long int) (fval);
376  if(identifier == "C"){
377  if(value > 0){
378  numParts=value;
379  } else {
380  throw "Invalid argument at " + tmp;
381  }
382  }
383  else if(identifier == "PF"){
384  stringstream stream(stringstream::in | stringstream::out);
385  stream << tmp;
386  getline(stream, paramFile, '=');
387  stream >> paramFile;
388  }
389  else {
390  throw "Invalid argument at " + tmp;
391  }
392 
393  }
394 
395 }
396 
397 void readGeoGenParams(string paramFileName, Teuchos::ParameterList &geoparams, const RCP<const Teuchos::Comm<int> > & comm){
398  std::string input = "";
399  char inp[25000];
400  for(int i = 0; i < 25000; ++i){
401  inp[i] = 0;
402  }
403  if(comm->getRank() == 0){
404  fstream inParam(paramFileName.c_str());
405 
406  std::string tmp = "";
407  getline (inParam,tmp);
408  while (!inParam.eof()){
409  if(tmp != ""){
410  tmp = trim_copy(tmp);
411  if(tmp != ""){
412  input += tmp + "\n";
413  }
414  }
415  getline (inParam,tmp);
416  }
417  inParam.close();
418  for (size_t i = 0; i < input.size(); ++i){
419  inp[i] = input[i];
420  }
421  }
422 
423 
424  int size = input.size();
425 
426  //MPI_Bcast(&size,1,MPI_INT, 0, MPI_COMM_WORLD);
427  //MPI_Bcast(inp,size,MPI_CHAR, 0, MPI_COMM_WORLD);
428  //Teuchos::broadcast<int, char>(comm, 0,inp);
429 
430  comm->broadcast(0, sizeof(int), (char*) &size);
431  comm->broadcast(0, size, inp);
432  //Teuchos::broadcast<int,string>(comm,0, &input);
433  istringstream inParam(inp);
434  string str;
435  getline (inParam,str);
436  while (!inParam.eof()){
437  if(str[0] != param_comment){
438  size_t pos = str.find('=');
439  if(pos == string::npos){
440  throw "Invalid Line:" + str + " in parameter file";
441  }
442  string paramname = trim_copy(str.substr(0,pos));
443  string paramvalue = trim_copy(str.substr(pos + 1));
444  geoparams.set(paramname, paramvalue);
445  }
446  getline (inParam,str);
447  }
448 }
449 
450 int main(int argc, char *argv[])
451 {
452  // MEMORY_CHECK(true, "Before initializing MPI");
453 
454  Teuchos::GlobalMPISession session(&argc, &argv, NULL);
455  RCP<const Comm<int> > comm = Teuchos::DefaultComm<int>::getComm();
456  int rank = comm->getRank();
457  int nprocs = comm->getSize();
458 
459  MEMORY_CHECK(rank==0 || rank==nprocs-1, "After initializing MPI");
460 
461  if (rank==0)
462  cout << "Number of processes: " << nprocs << endl;
463 
464  // Default values
465  //double numGlobalCoords = 1000;
466  //int nWeights = 0;
467  int debugLevel=2; // for timing
468  string memoryOn("memoryOn");
469  string memoryOff("memoryOff");
470  bool doMemory=false;
471 
472  int numGlobalParts = 512 * 512;
473  string paramFile = "p2.txt";
474  //getArgVals(argc, argv, numGlobalParts, paramFile);
475 
476  int dummyTimer=0;
477 
478 
479 
480  Teuchos::ParameterList geoparams("geo params");
481 
482  readGeoGenParams(paramFile, geoparams, comm);
483 #ifdef HAVE_ZOLTAN2_OMP
484  double begin = omp_get_wtime();
485 #endif
486  GeometricGenerator<zscalar_t, zlno_t, zgno_t, znode_t> *gg = new GeometricGenerator<zscalar_t, zlno_t, zgno_t, znode_t>(geoparams,comm);
487 #ifdef HAVE_ZOLTAN2_OMP
488  double end = omp_get_wtime();
489  cout << "GeometricGen Time:" << end - begin << endl;
490 #endif
491  int coord_dim = gg->getCoordinateDimension();
492  int nWeights = gg->getNumWeights();
493  zlno_t numLocalPoints = gg->getNumLocalCoords();
494  zgno_t numGlobalPoints = gg->getNumGlobalCoords();
495  zscalar_t **coords = new zscalar_t * [coord_dim];
496  for(int i = 0; i < coord_dim; ++i){
497  coords[i] = new zscalar_t[numLocalPoints];
498  }
499  gg->getLocalCoordinatesCopy(coords);
500  zscalar_t **weight = NULL;
501  if(nWeights){
502  weight= new zscalar_t * [nWeights];
503  for(int i = 0; i < nWeights; ++i){
504  weight[i] = new zscalar_t[numLocalPoints];
505  }
506  gg->getLocalWeightsCopy(weight);
507  }
508 
509  zgno_t globalSize = static_cast<zgno_t>(numGlobalPoints);
510  delete gg;
511 
512  cout << "coord_dim:" << coord_dim << " nWeights:" << nWeights << " numLocalPoints:" << numLocalPoints << " numGlobalPoints:" << numGlobalPoints << endl;
513 
514  CommandLineProcessor commandLine(false, true);
515  commandLine.setOption("size", &numGlobalPoints,
516  "Approximate number of global coordinates.");
517  commandLine.setOption("numParts", &numGlobalParts,
518  "Number of parts (default is one per proc).");
519  commandLine.setOption("numWeights", &nWeights,
520  "Number of weights per coordinate, zero implies uniform weights.");
521  commandLine.setOption("debug", &debugLevel, "Zoltan1 debug level");
522  commandLine.setOption("timers", &dummyTimer, "ignored");
523  commandLine.setOption(memoryOn.c_str(), memoryOff.c_str(), &doMemory,
524  "do memory profiling");
525 
526  string balanceCount("balance_object_count");
527  string balanceWeight("balance_object_weight");
528  string mcnorm1("multicriteria_minimize_total_weight");
529  string mcnorm2("multicriteria_balance_total_maximum");
530  string mcnorm3("multicriteria_minimize_maximum_weight");
531 
532  string objective(balanceWeight); // default
533 
534  string doc(balanceCount);
535  doc.append(": ignore weights\n");
536 
537  doc.append(balanceWeight);
538  doc.append(": balance on first weight\n");
539 
540  doc.append(mcnorm1);
541  doc.append(": given multiple weights, balance their total.\n");
542 
543  doc.append(mcnorm3);
544  doc.append(": given multiple weights, balance the maximum for each coordinate.\n");
545 
546  doc.append(mcnorm2);
547  doc.append(": given multiple weights, balance the L2 norm of the weights.\n");
548 
549  commandLine.setOption("objective", &objective, doc.c_str());
550 
551  CommandLineProcessor::EParseCommandLineReturn rc =
552  commandLine.parse(argc, argv);
553 
554  if (rc != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL){
555  if (rc == Teuchos::CommandLineProcessor::PARSE_HELP_PRINTED){
556  if (rank==0)
557  cout << "PASS" << endl;
558  return 1;
559  }
560  else{
561  if (rank==0)
562  cout << "FAIL" << endl;
563  return 0;
564  }
565  }
566 
567  //MEMORY_CHECK(doMemory && rank==0, "After processing parameters");
568 
569 
570 
571 
572  RCP<Tpetra::Map<zlno_t, zgno_t, znode_t> > mp = rcp(
573  new Tpetra::Map<zlno_t, zgno_t, znode_t> (numGlobalPoints, numLocalPoints, 0, comm));
574 
575  Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(coord_dim);
576  for (int i=0; i < coord_dim; i++){
577  if(numLocalPoints > 0){
578  Teuchos::ArrayView<const zscalar_t> a(coords[i], numLocalPoints);
579  coordView[i] = a;
580  } else{
581  Teuchos::ArrayView<const zscalar_t> a;
582  coordView[i] = a;
583  }
584  }
585 
586  coordinates = RCP< Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> >(
587  new Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>( mp, coordView.view(0, coord_dim), coord_dim));
588 
589  //typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> tMVector_t;
590 
591  RCP<const tMVector_t> coordsConst = Teuchos::rcp_const_cast<const tMVector_t>(coordinates);
592 
593 
594 
595  //coordinates = getMeshCoordinates(comm, globalSize);
596  size_t numLocalCoords = coordinates->getLocalLength();
597 
598 #if 0
599  comm->barrier();
600  for (int p=0; p < nprocs; p++){
601  if (p==rank){
602  cout << "Rank " << rank << ", " << numLocalCoords << "coords" << endl;
603  const zscalar_t *x = coordinates->getData(0).getRawPtr();
604  const zscalar_t *y = coordinates->getData(1).getRawPtr();
605  const zscalar_t *z = coordinates->getData(2).getRawPtr();
606  for (zlno_t i=0; i < numLocalCoords; i++)
607  cout << " " << x[i] << " " << y[i] << " " << z[i] << endl;
608  }
609  cout.flush();
610  comm->barrier();
611  }
612 #endif
613 
614  if (nWeights > 0){
615 
616  weights = arcp(new ArrayRCP<zscalar_t> [nWeights],
617  0, nWeights, true);
618  for(int i = 0; i < nWeights; ++i){
619  //weights[i] = ArrayRCP<zscalar_t>(weight[i]);
620  }
621  }
622 
623  MEMORY_CHECK(doMemory && rank==0, "After creating input");
624 
625  // Now call Zoltan to partition the problem.
626 
627  float ver;
628  int aok = Zoltan_Initialize(argc, argv, &ver);
629 
630  if (aok != 0){
631  printf("sorry...\n");
632  exit(0);
633  }
634 
635  struct Zoltan_Struct *zz;
636  zz = Zoltan_Create(MPI_COMM_WORLD);
637 
638  Zoltan_Set_Param(zz, "LB_METHOD", "RCB");
639  Zoltan_Set_Param(zz, "LB_APPROACH", "PARTITION");
640  Zoltan_Set_Param(zz, "CHECK_GEOM", "0");
641  Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1"); // compiled with ULONG option
642  Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "0");
643  Zoltan_Set_Param(zz, "RETURN_LISTS", "ALL");
644 
645  //Zoltan_Set_Param(zz, "FINAL_OUTPUT", "1");
646  Zoltan_Set_Param(zz, "IMBALANCE_TOL", "1.0");
647  std::ostringstream oss;
648  oss << numGlobalParts;
649  Zoltan_Set_Param(zz, "NUM_GLOBAL_PARTS", oss.str().c_str());
650  oss.str("");
651  oss << debugLevel;
652  Zoltan_Set_Param(zz, "DEBUG_LEVEL", oss.str().c_str());
653 
654  if (objective != balanceCount){
655  oss.str("");
656  oss << nWeights;
657  Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", oss.str().c_str());
658 
659  if (objective == mcnorm1)
660  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "1");
661  else if (objective == mcnorm2)
662  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "2");
663  else if (objective == mcnorm3)
664  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "3");
665  }
666  else{
667  Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", "0");
668  }
669 
670  Zoltan_Set_Num_Obj_Fn(zz, getNumObj, NULL);
671  Zoltan_Set_Obj_List_Fn(zz, getObjList,NULL);
672  Zoltan_Set_Num_Geom_Fn(zz, getDim, NULL);
673  Zoltan_Set_Geom_Multi_Fn(zz, getCoords, NULL);
674 
675  int changes, numGidEntries, numLidEntries, numImport, numExport;
676  zgno_t * importGlobalGids, * importLocalGids;
677  zgno_t * exportGlobalGids, * exportLocalGids;
678  int *importProcs, *importToPart, *exportProcs, *exportToPart;
679 
680  MEMORY_CHECK(doMemory && rank==0, "Before Zoltan_LB_Partition");
681 
682 
683  aok = Zoltan_LB_Partition(zz, /* input (all remaining fields are output) */
684  &changes, /* 1 if partitioning was changed, 0 otherwise */
685  &numGidEntries, /* Number of integers used for a global ID */
686  &numLidEntries, /* Number of integers used for a local ID */
687  &numImport, /* Number of vertices to be sent to me */
688  &importGlobalGids, /* Global IDs of vertices to be sent to me */
689  &importLocalGids, /* Local IDs of vertices to be sent to me */
690  &importProcs, /* Process rank for source of each incoming vertex */
691  &importToPart, /* New partition for each incoming vertex */
692  &numExport, /* Number of vertices I must send to other processes*/
693  &exportGlobalGids, /* Global IDs of the vertices I must send */
694  &exportLocalGids, /* Local IDs of the vertices I must send */
695  &exportProcs, /* Process to which I send each of the vertices */
696  &exportToPart); /* Partition to which each vertex will belong */
697 
698  MEMORY_CHECK(doMemory && rank==0, "After Zoltan_LB_Partition");
699  Zoltan_LB_Free_Part(importGlobalGids, importLocalGids, importProcs, importToPart);
700  Zoltan_LB_Free_Part(exportGlobalGids, exportLocalGids, exportProcs, exportToPart);
701  Zoltan_Destroy(&zz);
702  MEMORY_CHECK(doMemory && rank==0, "After Zoltan_Destroy");
703 
704  if (rank==0){
705  if (aok != 0)
706  std::cout << "FAIL" << std::endl;
707  else
708  std::cout << "PASS" << std::endl;
709  }
710 
711  return 0;
712 }
713 
714 #else
715 #include <iostream>
716 int main(int argc, char *argv[])
717 {
718  std::cout << "Test did not run due to faulty configuration." << std::endl;
719  std::cout << "FAIL" << std::endl;
720 }
721 #endif
string trim_right_copy(const string &s, const string &delimiters=" \f\n\r\t\v")
void readGeoGenParams(string paramFileName, Teuchos::ParameterList &geoparams, const RCP< const Teuchos::Comm< int > > &comm)
double zscalar_t
int zlno_t
Defines the PartitioningSolution class.
common code used by tests
int main(int argc, char *argv[])
Defines the XpetraMultiVectorAdapter.
string trim_copy(const string &s, const string &delimiters=" \f\n\r\t\v")
const char param_comment
Defines the PartitioningSolutionQuality class.
int zgno_t
string convert_to_string(char *args)
void getArgVals(int argc, char **argv, int &numParts, float &imbalance, string &pqParts, int &opt, std::string &fname, std::string &pfname, int &k, int &migration_check_option, int &migration_all_to_all_type, zscalar_t &migration_imbalance_cut_off, int &migration_processor_assignment_type, int &migration_doMigration_type, bool &test_boxes, bool &rectilinear)
string trim_left_copy(const string &s, const string &delimiters=" \f\n\r\t\v")
Tpetra::MultiVector< double, int, int > tMVector_t
bool getArgumentValue(string &argumentid, double &argumentValue, string argumentline)
Defines the PartitioningProblem class.
#define MEMORY_CHECK(iPrint, msg)