Zoltan2
scaling/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 
52 #include "Zoltan2_config.h"
53 #include <zoltan.h>
54 
55 #include <Zoltan2_Util.hpp>
56 #define MEMORY_CHECK(iPrint, msg) \
57  if (iPrint){ \
58  long kb = Zoltan2::getProcessKilobytes(); \
59  std::cout.width(10); \
60  std::cout.fill('*'); \
61  std::cout << kb << " KB, " << msg << std::endl; \
62  std::cout.width(0); \
63  std::cout.fill(' '); \
64  }
65 
66 #include <Teuchos_RCP.hpp>
67 #include <Teuchos_ArrayView.hpp>
68 #include <Teuchos_ParameterList.hpp>
69 #include <Teuchos_DefaultComm.hpp>
70 #include <Teuchos_Comm.hpp>
71 #include <Teuchos_CommHelpers.hpp>
72 
73 #include <Zoltan2_TestHelpers.hpp>
74 #include <Tpetra_MultiVector.hpp>
75 #include <Kokkos_DefaultNode.hpp>
76 #include <GeometricGenerator.hpp>
77 
78 #include <vector>
79 #include <string>
80 #include <ostream>
81 #include <sstream>
82 #include <fstream>
83 using std::string;
84 using std::vector;
85 using std::bad_alloc;
86 using Teuchos::RCP;
87 using Teuchos::rcp;
88 using Teuchos::Comm;
89 using Teuchos::ArrayView;
90 using Teuchos::CommandLineProcessor;
91 
92 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> tMVector_t;
93 typedef Tpetra::Map<zlno_t, zgno_t, znode_t> tMap_t;
94 typedef tMap_t::node_type znode_t;
95 
97 // Data structure for data
98 typedef struct dots {
99  vector<vector<float> > weights;
102 
103 const char param_comment = '#';
104 
106  const string& s,
107  const string& delimiters = " \f\n\r\t\v" )
108 {
109  return s.substr( 0, s.find_last_not_of( delimiters ) + 1 );
110 }
111 
113  const string& s,
114  const string& delimiters = " \f\n\r\t\v" )
115 {
116  return s.substr( s.find_first_not_of( delimiters ) );
117 }
118 
119 string trim_copy(
120  const string& s,
121  const string& delimiters = " \f\n\r\t\v" )
122 {
123  return trim_left_copy( trim_right_copy( s, delimiters ), delimiters );
124 }
125 
126 void readGeoGenParams(string paramFileName, Teuchos::ParameterList &geoparams, const RCP<const Teuchos::Comm<int> > & comm){
127  std::string input = "";
128  char inp[25000];
129  for(int i = 0; i < 25000; ++i){
130  inp[i] = 0;
131  }
132 
133  bool fail = false;
134  if(comm->getRank() == 0){
135 
136  std::fstream inParam(paramFileName.c_str());
137  if (inParam.fail())
138  {
139  fail = true;
140  }
141  if(!fail)
142  {
143  std::string tmp = "";
144  getline (inParam,tmp);
145  while (!inParam.eof()){
146  if(tmp != ""){
147  tmp = trim_copy(tmp);
148  if(tmp != ""){
149  input += tmp + "\n";
150  }
151  }
152  getline (inParam,tmp);
153  }
154  inParam.close();
155  for (size_t i = 0; i < input.size(); ++i){
156  inp[i] = input[i];
157  }
158  }
159  }
160 
161 
162 
163  int size = input.size();
164  if(fail){
165  size = -1;
166  }
167  comm->broadcast(0, sizeof(int), (char*) &size);
168  if(size == -1){
169  throw "File " + paramFileName + " cannot be opened.";
170  }
171  comm->broadcast(0, size, inp);
172  std::istringstream inParam(inp);
173  string str;
174  getline (inParam,str);
175  while (!inParam.eof()){
176  if(str[0] != param_comment){
177  size_t pos = str.find('=');
178  if(pos == string::npos){
179  throw "Invalid Line:" + str + " in parameter file";
180  }
181  string paramname = trim_copy(str.substr(0,pos));
182  string paramvalue = trim_copy(str.substr(pos + 1));
183  geoparams.set(paramname, paramvalue);
184  }
185  getline (inParam,str);
186  }
187 }
188 
189 
191 // Zoltan1 query functions
192 
193 int getNumObj(void *data, int *ierr)
194 {
195  *ierr = 0;
196  DOTS *dots = (DOTS *) data;
197  return dots->coordinates->getLocalLength();
198 }
199 
201 int getDim(void *data, int *ierr)
202 {
203  *ierr = 0;
204  DOTS *dots = (DOTS *) data;
205  int dim = dots->coordinates->getNumVectors();
206 
207  return dim;
208 }
209 
211 void getObjList(void *data, int numGid, int numLid,
212  ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids,
213  int num_wgts, float *obj_wgts, int *ierr)
214 {
215  *ierr = 0;
216  DOTS *dots = (DOTS *) data;
217 
218  size_t localLen = dots->coordinates->getLocalLength();
219  const zgno_t *ids =
220  dots->coordinates->getMap()->getNodeElementList().getRawPtr();
221 
222  if (sizeof(ZOLTAN_ID_TYPE) == sizeof(zgno_t))
223  memcpy(gids, ids, sizeof(ZOLTAN_ID_TYPE) * localLen);
224  else
225  for (size_t i=0; i < localLen; i++)
226  gids[i] = static_cast<ZOLTAN_ID_TYPE>(ids[i]);
227 
228  if (num_wgts > 0){
229  float *wgts = obj_wgts;
230  for (size_t i=0; i < localLen; i++)
231  for (int w=0; w < num_wgts; w++)
232  *wgts++ = dots->weights[w][i];
233  }
234 }
235 
237 void getCoords(void *data, int numGid, int numLid,
238  int numObj, ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids,
239  int dim, double *coords, int *ierr)
240 {
241  // I know that Zoltan asks for coordinates in gid order.
242  if (dim == 3){
243  *ierr = 0;
244  DOTS *dots = (DOTS *) data;
245  double *val = coords;
246  const zscalar_t *x = dots->coordinates->getData(0).getRawPtr();
247  const zscalar_t *y = dots->coordinates->getData(1).getRawPtr();
248  const zscalar_t *z = dots->coordinates->getData(2).getRawPtr();
249  for (int i=0; i < numObj; i++){
250  *val++ = static_cast<double>(x[i]);
251  *val++ = static_cast<double>(y[i]);
252  *val++ = static_cast<double>(z[i]);
253  }
254  }
255  else {
256  *ierr = 0;
257  DOTS *dots = (DOTS *) data;
258  double *val = coords;
259  const zscalar_t *x = dots->coordinates->getData(0).getRawPtr();
260  const zscalar_t *y = dots->coordinates->getData(1).getRawPtr();
261  for (int i=0; i < numObj; i++){
262  *val++ = static_cast<double>(x[i]);
263  *val++ = static_cast<double>(y[i]);
264  }
265 
266 
267  }
268 }
269 
270 
272 
278 };
279 
281  const RCP<const Teuchos::Comm<int> > & comm,
282  vector<float> &wgts, weightTypes how, float scale, int rank)
283 {
284  zlno_t len = wgts.size();
285  if (how == upDown){
286  float val = scale + rank%2;
287  for (zlno_t i=0; i < len; i++)
288  wgts[i] = val;
289  }
290  else if (how == roundRobin){
291  for (int i=0; i < 10; i++){
292  float val = (i + 10)*scale;
293  for (zlno_t j=i; j < len; j += 10)
294  wgts[j] = val;
295  }
296  }
297  else if (how == increasing){
298  float val = scale + rank;
299  for (zlno_t i=0; i < len; i++)
300  wgts[i] = val;
301  }
302 }
303 
305 /* Create a mesh of approximately the desired size.
306  *
307  * We want 3 dimensions close to equal in length.
308  */
310  const RCP<const Teuchos::Comm<int> > & comm,
311  zgno_t numGlobalCoords)
312 {
313  int rank = comm->getRank();
314  int nprocs = comm->getSize();
315 
316  double k = log(numGlobalCoords) / 3;
317  double xdimf = exp(k) + 0.5;
318  ssize_t xdim = static_cast<ssize_t>(floor(xdimf));
319  ssize_t ydim = xdim;
320  ssize_t zdim = numGlobalCoords / (xdim*ydim);
321  ssize_t num=xdim*ydim*zdim;
322  ssize_t diff = numGlobalCoords - num;
323  ssize_t newdiff = 0;
324 
325  while (diff > 0){
326  if (zdim > xdim && zdim > ydim){
327  zdim++;
328  newdiff = diff - (xdim*ydim);
329  if (newdiff < 0)
330  if (diff < -newdiff)
331  zdim--;
332  }
333  else if (ydim > xdim && ydim > zdim){
334  ydim++;
335  newdiff = diff - (xdim*zdim);
336  if (newdiff < 0)
337  if (diff < -newdiff)
338  ydim--;
339  }
340  else{
341  xdim++;
342  newdiff = diff - (ydim*zdim);
343  if (newdiff < 0)
344  if (diff < -newdiff)
345  xdim--;
346  }
347 
348  diff = newdiff;
349  }
350 
351  num=xdim*ydim*zdim;
352  diff = numGlobalCoords - num;
353  if (diff < 0)
354  diff /= -numGlobalCoords;
355  else
356  diff /= numGlobalCoords;
357 
358  if (rank == 0){
359  if (diff > .01)
360  std::cout << "Warning: Difference " << diff*100 << " percent" << std::endl;
361  std::cout << "Mesh size: " << xdim << "x" << ydim << "x" <<
362  zdim << ", " << num << " vertices." << std::endl;
363  }
364 
365  // Divide coordinates.
366 
367  ssize_t numLocalCoords = num / nprocs;
368  ssize_t leftOver = num % nprocs;
369  ssize_t gid0 = 0;
370 
371  if (rank <= leftOver)
372  gid0 = rank * (numLocalCoords+1);
373  else
374  gid0 = (leftOver * (numLocalCoords+1)) +
375  ((rank - leftOver) * numLocalCoords);
376 
377  if (rank < leftOver)
378  numLocalCoords++;
379 
380  ssize_t gid1 = gid0 + numLocalCoords;
381 
382  zgno_t *ids = new zgno_t[numLocalCoords];
383  if (!ids)
384  throw bad_alloc();
385  ArrayView<zgno_t> idArray(ids, numLocalCoords);
386  zgno_t *idptr = ids;
387 
388  for (ssize_t i=gid0; i < gid1; i++)
389  *idptr++ = zgno_t(i);
390 
391  RCP<const tMap_t> idMap = rcp(new tMap_t(num, idArray, 0, comm));
392 
393  delete [] ids;
394 
395  // Create a Tpetra::MultiVector of coordinates.
396 
397  zscalar_t *x = new zscalar_t [numLocalCoords*3];
398  if (!x) throw bad_alloc();
399 
400  zscalar_t *y = x + numLocalCoords;
401  zscalar_t *z = y + numLocalCoords;
402 
403  zgno_t xStart = 0;
404  zgno_t yStart = 0;
405  zgno_t xyPlane = xdim*ydim;
406  zgno_t zStart = gid0 / xyPlane;
407  zgno_t rem = gid0 % xyPlane;
408  if (rem > 0){
409  yStart = rem / xdim;
410  xStart = rem % xdim;
411  }
412 
413  zlno_t next = 0;
414  for (zscalar_t zval=zStart; next < numLocalCoords && zval < zdim; zval+=1.){
415  for (zscalar_t yval=yStart; next < numLocalCoords && yval < ydim; yval+=1.){
416  for (zscalar_t xval=xStart; next < numLocalCoords && xval < xdim;xval+=1.){
417  x[next] = xval;
418  y[next] = yval;
419  z[next] = zval;
420  next++;
421  }
422  xStart = 0;
423  }
424  yStart = 0;
425  }
426 
427  ArrayView<const zscalar_t> xArray(x, numLocalCoords*3);
428  tMVector_t *dots = new tMVector_t(idMap, xArray, numLocalCoords, 3);
429 
430  delete [] x;
431  return dots;
432 }
433 
434 
435 int main(int narg, char *arg[])
436 {
437  // MEMORY_CHECK(true, "Before initializing MPI");
438 
439  Tpetra::ScopeGuard tscope(&narg, &arg);
440  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
441  int rank = comm->getRank();
442  int nprocs = comm->getSize();
443  DOTS dots;
444 
445  MEMORY_CHECK(rank==0 || rank==nprocs-1, "After initializing MPI");
446 
447  if (rank==0)
448  std::cout << "Number of processes: " << nprocs << std::endl;
449 
450  // Default values
451  zgno_t numGlobalCoords = 1000;
452  int nWeights = 0;
453  int debugLevel=2;
454  string memoryOn("memoryOn");
455  string memoryOff("memoryOff");
456  bool doMemory=false;
457  int numGlobalParts = nprocs;
458  int dummyTimer=0;
459  bool remap=0;
460 
461  string balanceCount("balance_object_count");
462  string balanceWeight("balance_object_weight");
463  string mcnorm1("multicriteria_minimize_total_weight");
464  string mcnorm2("multicriteria_balance_total_maximum");
465  string mcnorm3("multicriteria_minimize_maximum_weight");
466  string objective(balanceWeight); // default
467 
468  // Process command line input
469  CommandLineProcessor commandLine(false, true);
470  //commandLine.setOption("size", &numGlobalCoords,
471  // "Approximate number of global coordinates.");
472  int input_option = 0;
473  commandLine.setOption("input_option", &input_option,
474  "whether to use mesh creation, geometric generator, or file input");
475  string inputFile = "";
476 
477  commandLine.setOption("input_file", &inputFile,
478  "the input file for geometric generator or file input");
479 
480 
481  commandLine.setOption("size", &numGlobalCoords,
482  "Approximate number of global coordinates.");
483  commandLine.setOption("numParts", &numGlobalParts,
484  "Number of parts (default is one per proc).");
485  commandLine.setOption("nWeights", &nWeights,
486  "Number of weights per coordinate, zero implies uniform weights.");
487  commandLine.setOption("debug", &debugLevel, "Zoltan1 debug level");
488  commandLine.setOption("remap", "no-remap", &remap,
489  "Zoltan1 REMAP parameter; disabled by default for scalability testing");
490  commandLine.setOption("timers", &dummyTimer, "ignored");
491  commandLine.setOption(memoryOn.c_str(), memoryOff.c_str(), &doMemory,
492  "do memory profiling");
493 
494  string doc(balanceCount);
495  doc.append(": ignore weights\n");
496  doc.append(balanceWeight);
497  doc.append(": balance on first weight\n");
498  doc.append(mcnorm1);
499  doc.append(": given multiple weights, balance their total.\n");
500  doc.append(mcnorm3);
501  doc.append(": given multiple weights, "
502  "balance the maximum for each coordinate.\n");
503  doc.append(mcnorm2);
504  doc.append(": given multiple weights, balance the L2 norm of the weights.\n");
505  commandLine.setOption("objective", &objective, doc.c_str());
506 
507  CommandLineProcessor::EParseCommandLineReturn rc =
508  commandLine.parse(narg, arg);
509 
510 
511 
512  if (rc != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
513  if (rc == Teuchos::CommandLineProcessor::PARSE_HELP_PRINTED) {
514  if (rank==0) std::cout << "PASS" << std::endl;
515  return 1;
516  }
517  else {
518  if (rank==0) std::cout << "FAIL" << std::endl;
519  return 0;
520  }
521  }
522 
523  //MEMORY_CHECK(doMemory && rank==0, "After processing parameters");
524 
525  // Create the data structure
526  size_t numLocalCoords = 0;
527  if (input_option == 0){
528  dots.coordinates = makeMeshCoordinates(comm, numGlobalCoords);
529  numLocalCoords = dots.coordinates->getLocalLength();
530 
531 #if 0
532  comm->barrier();
533  for (int p=0; p < nprocs; p++){
534  if (p==rank){
535  std::cout << "Rank " << rank << ", " << numLocalCoords << "coords" << std::endl;
536  const zscalar_t *x = coordinates->getData(0).getRawPtr();
537  const zscalar_t *y = coordinates->getData(1).getRawPtr();
538  const zscalar_t *z = coordinates->getData(2).getRawPtr();
539  for (zlno_t i=0; i < numLocalCoords; i++)
540  std::cout << " " << x[i] << " " << y[i] << " " << z[i] << std::endl;
541  }
542  std::cout.flush();
543  comm->barrier();
544  }
545 #endif
546 
547  if (nWeights > 0){
548 
549  dots.weights.resize(nWeights);
550 
551  int wt = 0;
552  float scale = 1.0;
553  for (int i=0; i < nWeights; i++){
554  dots.weights[i].resize(numLocalCoords);
555  makeWeights(comm, dots.weights[i], weightTypes(wt++), scale, rank);
556 
557  if (wt == numWeightTypes){
558  wt = 0;
559  scale++;
560  }
561  }
562  }
563  }
564  else if(input_option == 1){
565  Teuchos::ParameterList geoparams("geo params");
566  readGeoGenParams(inputFile, geoparams, comm);
568 
569  int coord_dim = gg->getCoordinateDimension();
570  nWeights = gg->getNumWeights();
571  numLocalCoords = gg->getNumLocalCoords();
572  numGlobalCoords = gg->getNumGlobalCoords();
573  zscalar_t **coords = new zscalar_t * [coord_dim];
574  for(int i = 0; i < coord_dim; ++i){
575  coords[i] = new zscalar_t[numLocalCoords];
576  }
577  gg->getLocalCoordinatesCopy(coords);
578  zscalar_t **weight = NULL;
579  if(nWeights){
580  weight= new zscalar_t * [nWeights];
581  for(int i = 0; i < nWeights; ++i){
582  weight[i] = new zscalar_t[numLocalCoords];
583  }
584  gg->getLocalWeightsCopy(weight);
585  }
586 
587  delete gg;
588 
589  RCP<Tpetra::Map<zlno_t, zgno_t, znode_t> > mp = rcp(
590  new Tpetra::Map<zlno_t, zgno_t, znode_t> (numGlobalCoords, numLocalCoords, 0, comm));
591 
592  Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(coord_dim);
593  for (int i=0; i < coord_dim; i++){
594  if(numLocalCoords > 0){
595  Teuchos::ArrayView<const zscalar_t> a(coords[i], numLocalCoords);
596  coordView[i] = a;
597  } else{
598  Teuchos::ArrayView<const zscalar_t> a;
599  coordView[i] = a;
600  }
601  }
602 
603  tMVector_t *tmVector = new tMVector_t( mp, coordView.view(0, coord_dim), coord_dim);
604 
605  dots.coordinates = tmVector;
606  dots.weights.resize(nWeights);
607 
608  if(nWeights){
609  for (int i = 0; i < nWeights;++i){
610  for (zlno_t j = 0; j < zlno_t(numLocalCoords); ++j){
611  dots.weights[i].push_back(weight[i][j]);
612  }
613  }
614  }
615  if(nWeights){
616  for(int i = 0; i < nWeights; ++i)
617  delete [] weight[i];
618  delete [] weight;
619  }
620  }
621  else {
622 
623  UserInputForTests uinput(testDataFilePath, inputFile, comm, true);
624  RCP<tMVector_t> coords = uinput.getUICoordinates();
625  tMVector_t *newMulti = new tMVector_t(*coords);
626  dots.coordinates = newMulti;
627  numLocalCoords = coords->getLocalLength();
628  numGlobalCoords = coords->getGlobalLength();
629  }
630 
631  MEMORY_CHECK(doMemory && rank==0, "After creating input");
632 
633  // Now call Zoltan to partition the problem.
634 
635  float ver;
636  int aok = Zoltan_Initialize(narg, arg, &ver);
637 
638  if (aok != 0){
639  printf("Zoltan_Initialize failed\n");
640  exit(0);
641  }
642 
643  struct Zoltan_Struct *zz;
644  zz = Zoltan_Create(MPI_COMM_WORLD);
645 
646  Zoltan_Set_Param(zz, "LB_METHOD", "RCB");
647  Zoltan_Set_Param(zz, "LB_APPROACH", "PARTITION");
648  Zoltan_Set_Param(zz, "CHECK_GEOM", "0");
649  Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1");
650  Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "0");
651  Zoltan_Set_Param(zz, "RETURN_LISTS", "PART");
652  std::ostringstream oss;
653  oss << numGlobalParts;
654  Zoltan_Set_Param(zz, "NUM_GLOBAL_PARTS", oss.str().c_str());
655  oss.str("");
656  oss << debugLevel;
657  Zoltan_Set_Param(zz, "DEBUG_LEVEL", oss.str().c_str());
658 
659  if (remap)
660  Zoltan_Set_Param(zz, "REMAP", "1");
661  else
662  Zoltan_Set_Param(zz, "REMAP", "0");
663 
664  if (objective != balanceCount){
665  oss.str("");
666  oss << nWeights;
667  Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", oss.str().c_str());
668 
669  if (objective == mcnorm1)
670  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "1");
671  else if (objective == mcnorm2)
672  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "2");
673  else if (objective == mcnorm3)
674  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "3");
675  }
676  else{
677  Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", "0");
678  }
679 
680  Zoltan_Set_Num_Obj_Fn(zz, getNumObj, &dots);
681  Zoltan_Set_Obj_List_Fn(zz, getObjList, &dots);
682  Zoltan_Set_Num_Geom_Fn(zz, getDim, &dots);
683  Zoltan_Set_Geom_Multi_Fn(zz, getCoords, &dots);
684 
685  int changes, numGidEntries, numLidEntries, numImport, numExport;
686  ZOLTAN_ID_PTR importGlobalGids, importLocalGids;
687  ZOLTAN_ID_PTR exportGlobalGids, exportLocalGids;
688  int *importProcs, *importToPart, *exportProcs, *exportToPart;
689 
690  MEMORY_CHECK(doMemory && rank==0, "Before Zoltan_LB_Partition");
691 
692  if (rank == 0) std::cout << "Calling Zoltan_LB_Partition" << std::endl;
693  aok = Zoltan_LB_Partition(zz, &changes, &numGidEntries, &numLidEntries,
694  &numImport, &importGlobalGids, &importLocalGids,
695  &importProcs, &importToPart,
696  &numExport, &exportGlobalGids, &exportLocalGids,
697  &exportProcs, &exportToPart);
698  if (rank == 0) std::cout << "Returned from Zoltan_LB_Partition" << std::endl;
699 
700  MEMORY_CHECK(doMemory && rank==0, "After Zoltan_LB_Partition");
701 
702  /* Print the load-balance stats here */
703 
704  float *sumWgtPerPart = new float[numGlobalParts];
705  float *gsumWgtPerPart = new float[numGlobalParts];
706  for (int i = 0; i < numGlobalParts; i++) sumWgtPerPart[i] = 0.;
707 
708  for (size_t i = 0; i < numLocalCoords; i++)
709  sumWgtPerPart[exportToPart[i]] += (nWeights ? dots.weights[0][i]: 1.);
710 
711  Teuchos::reduceAll<int, float>(*comm, Teuchos::REDUCE_SUM, numGlobalParts,
712  sumWgtPerPart, gsumWgtPerPart);
713 
714  float maxSumWgtPerPart = 0.;
715  float minSumWgtPerPart = std::numeric_limits<float>::max();
716  float totWgt = 0.;
717  int maxSumWgtPart=0, minSumWgtPart=0;
718  for (int i = 0; i < numGlobalParts; i++) {
719  if (gsumWgtPerPart[i] > maxSumWgtPerPart) {
720  maxSumWgtPerPart = gsumWgtPerPart[i];
721  maxSumWgtPart = i;
722  }
723  if (gsumWgtPerPart[i] < minSumWgtPerPart) {
724  minSumWgtPerPart = gsumWgtPerPart[i];
725  minSumWgtPart = i;
726  }
727  totWgt += gsumWgtPerPart[i];
728  }
729 
730  if (rank == 0)
731  std::cout << std::endl << std::endl
732  << "Part loads (per part for " << numGlobalParts << " parts):"
733  << std::endl
734  << " min = " << minSumWgtPerPart
735  << " in part " << minSumWgtPart << std::endl
736  << " max = " << maxSumWgtPerPart
737  << " in part " << maxSumWgtPart << std::endl
738  << " tot = " << totWgt << std::endl
739  << " avg = " << totWgt / numGlobalParts
740  << std::endl << std::endl << std::endl;
741 
742  delete [] sumWgtPerPart;
743  delete [] gsumWgtPerPart;
744 
745  Zoltan_Destroy(&zz);
746  MEMORY_CHECK(doMemory && rank==0, "After Zoltan_Destroy");
747 
748  delete dots.coordinates;
749  for (int i = 0; i < nWeights; i++)
750  dots.weights[i].clear();
751  dots.weights.clear();
752 
753  MEMORY_CHECK(doMemory && rank==0, "After destroying input");
754 
755  if (rank==0){
756  if (aok != 0)
757  std::cout << "FAIL" << std::endl;
758  else
759  std::cout << "PASS" << std::endl;
760  }
761 
762  return 0;
763 }
common code used by tests
float zscalar_t
Tpetra::Map ::local_ordinal_type zlno_t
std::string testDataFilePath(".")
Tpetra::Map ::global_ordinal_type zgno_t
A gathering of useful namespace methods.
RCP< tMVector_t > getUICoordinates()
static const std::string fail
static RCP< tMVector_t > coordinates
tMap_t::node_type znode_t
#define MEMORY_CHECK(iPrint, msg)
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tMVector_t
void readGeoGenParams(string paramFileName, Teuchos::ParameterList &geoparams, const RCP< const Teuchos::Comm< int > > &comm)
string trim_left_copy(const string &s, const string &delimiters=" \f\n\r\t\v")
tMVector_t * makeMeshCoordinates(const RCP< const Teuchos::Comm< int > > &comm, zgno_t numGlobalCoords)
string trim_right_copy(const string &s, const string &delimiters=" \f\n\r\t\v")
string trim_copy(const string &s, const string &delimiters=" \f\n\r\t\v")
int main(int narg, char *arg[])
const char param_comment
struct dots DOTS
void makeWeights(const RCP< const Teuchos::Comm< int > > &comm, vector< float > &wgts, weightTypes how, float scale, int rank)
void getCoords(void *data, int numGid, int numLid, int numObj, ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids, int dim, double *coords, int *ierr)
void getObjList(void *data, int numGid, int numLid, ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids, int num_wgts, float *obj_wgts, int *ierr)
int getDim(void *data, int *ierr)
Tpetra::Map< zlno_t, zgno_t, znode_t > tMap_t
int getNumObj(void *data, int *ierr)
vector< vector< float > > weights
tMVector_t * coordinates
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tMVector_t