Zoltan2
GeometricGenerator.hpp
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 #ifndef GEOMETRICGENERATOR
46 #define GEOMETRICGENERATOR
47 
48 #include <Teuchos_Comm.hpp>
49 #include <Teuchos_ParameterList.hpp>
50 #include <Teuchos_FilteredIterator.hpp>
51 #include <Teuchos_ParameterEntry.hpp>
52 #include <iostream>
53 #include <ctime>
54 #include <limits>
55 #include <climits>
56 #include <string>
57 #include <cstdlib>
58 #include <sstream>
59 #include <fstream>
60 #include <Tpetra_MultiVector_decl.hpp>
63 #include <Teuchos_ArrayViewDecl.hpp>
64 #include <Teuchos_RCP.hpp>
65 #include <Tpetra_Distributor.hpp>
67 
68 
69 //#define HAVE_ZOLTAN2_ZOLTAN
70 #ifdef HAVE_ZOLTAN2_ZOLTAN
71 #include <zoltan.h>
72 #endif
73 
74 #ifdef _MSC_VER
75 #define NOMINMAX
76 #include <windows.h>
77 #endif
78 
79 using Teuchos::CommandLineProcessor;
80 using Teuchos::ArrayView;
81 using std::string;
82 using Teuchos::ArrayRCP;
83 
84 namespace GeometricGen{
85 #define CATCH_EXCEPTIONS(pp) \
86  catch (std::runtime_error &e) { \
87  cout << "Runtime exception returned from " << pp << ": " \
88  << e.what() << " FAIL" << endl; \
89  return -1; \
90  } \
91  catch (std::logic_error &e) { \
92  cout << "Logic exception returned from " << pp << ": " \
93  << e.what() << " FAIL" << endl; \
94  return -1; \
95  } \
96  catch (std::bad_alloc &e) { \
97  cout << "Bad_alloc exception returned from " << pp << ": " \
98  << e.what() << " FAIL" << endl; \
99  return -1; \
100  } \
101  catch (std::exception &e) { \
102  cout << "Unknown exception returned from " << pp << ": " \
103  << e.what() << " FAIL" << endl; \
104  return -1; \
105  }
106 
107 
108 
109 
110 #ifdef HAVE_ZOLTAN2_ZOLTAN
111 
112 
113 template <typename tMVector_t>
114 class DOTS{
115 public:
116  std::vector<std::vector<float> > weights;
117  tMVector_t *coordinates;
118 };
119 
120 template <typename tMVector_t>
121 int getNumObj(void *data, int *ierr)
122 {
123  *ierr = 0;
124  DOTS<tMVector_t> *dots_ = (DOTS<tMVector_t> *) data;
125  return dots_->coordinates->getLocalLength();
126 }
128 template <typename tMVector_t>
129 void getCoords(void *data, int numGid, int numLid,
130  int numObj, ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids,
131  int dim, double *coords_, int *ierr)
132 {
133  typedef typename tMVector_t::scalar_type scalar_t;
134 
135  // I know that Zoltan asks for coordinates in gid order.
136  if (dim == 3){
137  *ierr = 0;
138  DOTS<tMVector_t> *dots_ = (DOTS<tMVector_t> *) data;
139  double *val = coords_;
140  const scalar_t *x = dots_->coordinates->getData(0).getRawPtr();
141  const scalar_t *y = dots_->coordinates->getData(1).getRawPtr();
142  const scalar_t *z = dots_->coordinates->getData(2).getRawPtr();
143  for (int i=0; i < numObj; i++){
144  *val++ = static_cast<double>(x[i]);
145  *val++ = static_cast<double>(y[i]);
146  *val++ = static_cast<double>(z[i]);
147  }
148  }
149  else {
150  *ierr = 0;
151  DOTS<tMVector_t> *dots_ = (DOTS<tMVector_t> *) data;
152  double *val = coords_;
153  const scalar_t *x = dots_->coordinates->getData(0).getRawPtr();
154  const scalar_t *y = dots_->coordinates->getData(1).getRawPtr();
155  for (int i=0; i < numObj; i++){
156  *val++ = static_cast<double>(x[i]);
157  *val++ = static_cast<double>(y[i]);
158  }
159  }
160 }
161 
162 template <typename tMVector_t>
163 int getDim(void *data, int *ierr)
164 {
165  *ierr = 0;
166  DOTS<tMVector_t> *dots_ = (DOTS<tMVector_t> *) data;
167  int dim = dots_->coordinates->getNumVectors();
168 
169  return dim;
170 }
171 
173 template <typename tMVector_t>
174 void getObjList(void *data, int numGid, int numLid,
175  ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids,
176  int num_wgts, float *obj_wgts, int *ierr)
177 {
178  typedef typename tMVector_t::global_ordinal_type gno_t;
179 
180  *ierr = 0;
181  DOTS<tMVector_t> *dots_ = (DOTS<tMVector_t> *) data;
182 
183  size_t localLen = dots_->coordinates->getLocalLength();
184  const gno_t *ids =
185  dots_->coordinates->getMap()->getNodeElementList().getRawPtr();
186 
187  if (sizeof(ZOLTAN_ID_TYPE) == sizeof(gno_t))
188  memcpy(gids, ids, sizeof(ZOLTAN_ID_TYPE) * localLen);
189  else
190  for (size_t i=0; i < localLen; i++)
191  gids[i] = static_cast<ZOLTAN_ID_TYPE>(ids[i]);
192 
193  if (num_wgts > 0){
194  float *wgts = obj_wgts;
195  for (size_t i=0; i < localLen; i++)
196  for (int w=0; w < num_wgts; w++)
197  *wgts++ = dots_->weights[w][i];
198  }
199 }
200 #endif
201 
202 
204 const std::string shapes[] = {"SQUARE", "RECTANGLE", "CIRCLE", "CUBE", "RECTANGULAR_PRISM", "SPHERE"};
205 #define SHAPE_COUNT 6
206 
208 const std::string distribution[] = {"distribution", "uniform"};
209 #define DISTRIBUTION_COUNT 2
210 
211 #define HOLE_ALLOC_STEP 10
212 #define MAX_WEIGHT_DIM 10
213 #define INVALID(STR) "Invalid argument at " + STR
214 #define INVALIDSHAPE(STR, DIM) "Invalid shape name " + STR + " for " + DIM + ".\nValid shapes are \"SQUARE\", \"RECTANGLE\", \"CIRCLE\" for 2D, and \"CUBE\", \"RECTANGULAR_PRISM\", \"SPHERE\" for 3D"
215 
216 #define INVALID_SHAPE_ARG(SHAPE, REQUIRED) "Invalid argument count for shape " + SHAPE + ". Requires " + REQUIRED + " argument(s)."
217 #define MAX_ITER_ALLOWED 500
218 
219 const std::string weight_distribution_string = "WeightDistribution-";
220 
221 template <typename T>
223  T x;
224  T y;
225  T z;
226  /*
227  bool isInArea(int dim, T *minCoords, T *maxCoords){
228  dim = min(3, dim);
229  for (int i = 0; i < dim; ++i){
230  bool maybe = true;
231  switch(i){
232  case 0:
233  maybe = x < maxCoords[i] && x > maxCoords[i];
234  break;
235  case 1:
236  maybe = y < maxCoords[i] && y > maxCoords[i];
237  break;
238  case 2:
239  maybe = z < maxCoords[i] && z > maxCoords[i];
240  break;
241  }
242  if (!maybe){
243  return false;
244  }
245  }
246  return true;
247  }
248  */
250  x = 0;y=0;z=0;
251  }
252 };
253 
254 template <typename T>
255 class Hole{
256 public:
258  T edge1, edge2, edge3;
259  Hole(CoordinatePoint<T> center_, T edge1_, T edge2_, T edge3_){
260  this->center.x = center_.x;
261  this->center.y = center_.y;
262  this->center.z = center_.z;
263  this->edge1 = edge1_;
264  this->edge2 = edge2_;
265  this->edge3 = edge3_;
266  }
267  virtual bool isInArea(CoordinatePoint<T> dot) = 0;
268  virtual ~Hole(){}
269 };
270 
271 template <typename T>
272 class SquareHole:public Hole<T>{
273 public:
274  SquareHole(CoordinatePoint<T> center_ , T edge_): Hole<T>(center_, edge_, 0 , 0){}
275  virtual ~SquareHole(){}
276  virtual bool isInArea(CoordinatePoint<T> dot){
277  return fabs(dot.x - this->center.x) < this->edge1 / 2 && fabs(dot.y - this->center.y) < this->edge1 / 2;
278  }
279 };
280 
281 template <typename T>
282 class RectangleHole:public Hole<T>{
283 public:
284  RectangleHole(CoordinatePoint<T> center_ , T edge_x_, T edge_y_): Hole<T>(center_, edge_x_, edge_y_, 0){}
285  virtual bool isInArea(CoordinatePoint<T> dot){
286  return fabs(dot.x - this->center.x) < this->edge1 / 2 && fabs(dot.y - this->center.y) < this->edge2 / 2;
287  }
288  virtual ~RectangleHole(){}
289 };
290 
291 template <typename T>
292 class CircleHole:public Hole<T>{
293 public:
294  CircleHole(CoordinatePoint<T> center_ , T edge_): Hole<T>(center_, edge_, 0 , 0){}
295  virtual ~CircleHole(){}
296  virtual bool isInArea(CoordinatePoint<T> dot){
297  return (dot.x - this->center.x)*(dot.x - this->center.x) + (dot.y - this->center.y) * (dot.y - this->center.y) < this->edge1 * this->edge1;
298  }
299 };
300 
301 template <typename T>
302 class CubeHole:public Hole<T>{
303 public:
304  CubeHole(CoordinatePoint<T> center_ , T edge_): Hole<T>(center_, edge_, 0 , 0){}
305  virtual ~CubeHole(){}
306  virtual bool isInArea(CoordinatePoint<T> dot){
307  return fabs(dot.x - this->center.x) < this->edge1 / 2 && fabs(dot.y - this->center.y) < this->edge1 / 2 && fabs(dot.z - this->center.z) < this->edge1 / 2;
308  }
309 };
310 
311 template <typename T>
312 class RectangularPrismHole:public Hole<T>{
313 public:
314  RectangularPrismHole(CoordinatePoint<T> center_ , T edge_x_, T edge_y_, T edge_z_): Hole<T>(center_, edge_x_, edge_y_, edge_z_){}
316  virtual bool isInArea(CoordinatePoint<T> dot){
317  return fabs(dot.x - this->center.x) < this->edge1 / 2 && fabs(dot.y - this->center.y) < this->edge2 / 2 && fabs(dot.z - this->center.z) < this->edge3 / 2;
318  }
319 };
320 
321 template <typename T>
322 class SphereHole:public Hole<T>{
323 public:
324  SphereHole(CoordinatePoint<T> center_ , T edge_): Hole<T>(center_, edge_, 0 , 0){}
325  virtual ~SphereHole(){}
326  virtual bool isInArea(CoordinatePoint<T> dot){
327  return fabs((dot.x - this->center.x) * (dot.x - this->center.x) * (dot.x - this->center.x)) +
328  fabs((dot.y - this->center.y) * (dot.y - this->center.y) * (dot.y - this->center.y)) +
329  fabs((dot.z - this->center.z) * (dot.z - this->center.z) * (dot.z - this->center.z))
330  <
331  this->edge1 * this->edge1 * this->edge1;
332  }
333 };
334 
335 template <typename T, typename weighttype>
337 public:
338  virtual weighttype getWeight(CoordinatePoint<T> P)=0;
339  virtual weighttype get1DWeight(T x)=0;
340  virtual weighttype get2DWeight(T x, T y)=0;
341  virtual weighttype get3DWeight(T x, T y, T z)=0;
343  virtual ~WeightDistribution(){};
344 };
345 
364 template <typename T, typename weighttype>
365 class SteppedEquation:public WeightDistribution<T, weighttype>{
366  T a1,a2,a3;
367  T b1,b2,b3;
368  T c;
369  T x1,y1,z1;
370 
371  T *steps;
372  weighttype *values;
373  int stepCount;
374 public:
375  SteppedEquation(T a1_, T a2_, T a3_, T b1_, T b2_, T b3_, T c_, T x1_, T y1_, T z1_, T *steps_, T *values_, int stepCount_):WeightDistribution<T,weighttype>(){
376  this->a1 = a1_;
377  this->a2 = a2_;
378  this->a3 = a3_;
379  this->b1 = b1_;
380  this->b2 = b2_;
381  this->b3 = b3_;
382  this->c = c_;
383  this->x1 = x1_;
384  this->y1 = y1_;
385  this->z1 = z1_;
386 
387 
388  this->stepCount = stepCount_;
389  if(this->stepCount > 0){
390  this->steps = new T[this->stepCount];
391  this->values = new T[this->stepCount + 1];
392 
393  for (int i = 0; i < this->stepCount; ++i){
394  this->steps[i] = steps_[i];
395  this->values[i] = values_[i];
396  }
397  this->values[this->stepCount] = values_[this->stepCount];
398  }
399 
400  }
401 
402  virtual ~SteppedEquation(){
403  if(this->stepCount > 0){
404  delete [] this->steps;
405  delete [] this->values;
406  }
407  }
408 
409 
410  virtual weighttype get1DWeight(T x){
411  T expressionRes = this->a1 * pow( (x - this->x1), b1) + c;
412  if(this->stepCount > 0){
413  for (int i = 0; i < this->stepCount; ++i){
414  if (expressionRes < this->steps[i]) return this->values[i];
415  }
416  return this->values[this->stepCount];
417  }
418  else {
419  return weighttype(expressionRes);
420  }
421  }
422 
423  virtual weighttype get2DWeight(T x, T y){
424  T expressionRes = this->a1 * pow( (x - this->x1), b1) + this->a2 * pow( (y - this->y1), b2) + c;
425  if(this->stepCount > 0){
426  for (int i = 0; i < this->stepCount; ++i){
427  if (expressionRes < this->steps[i]) return this->values[i];
428  }
429  return this->values[this->stepCount];
430  }
431  else {
432  return weighttype(expressionRes);
433  }
434  }
435 
436  void print (T x, T y, T z){
437  cout << this->a1 << "*" << "math.pow( (" << x << "- " << this->x1 << " ), " << b1 <<")" << "+" << this->a2<< "*"<< "math.pow( (" << y << "-" << this->y1 << "), " <<
438  b2 << " ) + " << this->a3 << " * math.pow( (" << z << "-" << this->z1 << "), " << b3 << ")+ " << c << " == " <<
439  this->a1 * pow( (x - this->x1), b1) + this->a2 * pow( (y - this->y1), b2) + this->a3 * pow( (z - this->z1), b3) + c << endl;
440 
441  }
442 
443  virtual weighttype get3DWeight(T x, T y, T z){
444  T expressionRes = this->a1 * pow( (x - this->x1), b1) + this->a2 * pow( (y - this->y1), b2) + this->a3 * pow( (z - this->z1), b3) + this->c;
445 
446  //this->print(x,y,z);
447  if(this->stepCount > 0){
448  for (int i = 0; i < this->stepCount; ++i){
449  if (expressionRes < this->steps[i]) {
450  //cout << "0exp:" << expressionRes << " step:" << steps[i] << " value:" << values[i] << endl;
451  return this->values[i];
452  }
453  }
454 
455  //cout << "1exp:" << expressionRes << " step:" << steps[stepCount] << " value:" << values[stepCount] << endl;
456  return this->values[this->stepCount];
457  }
458  else {
459  return weighttype(expressionRes);
460  }
461  }
462 
463  virtual weighttype getWeight(CoordinatePoint<T> p){
464  T expressionRes = this->a1 * pow( (p.x - this->x1), b1) + this->a2 * pow( (p.y - this->y1), b2) + this->a3 * pow( (p.z - this->z1), b3);
465  if(this->stepCount > 0){
466  for (int i = 0; i < this->stepCount; ++i){
467  if (expressionRes < this->steps[i]) return this->values[i];
468  }
469  return this->values[this->stepCount];
470  }
471  else {
472  return weighttype(expressionRes);
473  }
474  }
475 };
476 
477 template <typename T, typename lno_t, typename gno_t>
479 public:
480  gno_t numPoints;
482  lno_t requested;
486 
487  CoordinateDistribution(gno_t np_, int dim, int wSize) :
488  numPoints(np_), dimension(dim), requested(0), assignedPrevious(0),
489  worldSize(wSize){}
490 
491  virtual CoordinatePoint<T> getPoint(gno_t point_index, unsigned int &state) = 0;
492  virtual T getXCenter() = 0;
493  virtual T getXRadius() =0;
494 
495  void GetPoints(lno_t requestedPointcount, CoordinatePoint<T> *points /*preallocated sized numPoints*/,
496  Hole<T> **holes, lno_t holeCount,
497  float *sharedRatios_, int myRank){
498 
499  for (int i = 0; i < myRank; ++i){
500  //cout << "me:" << myRank << " i:" << i << " s:" << sharedRatios_[i]<< endl;
501  this->assignedPrevious += lno_t(sharedRatios_[i] * this->numPoints);
502  if (i < this->numPoints % this->worldSize){
503  this->assignedPrevious += 1;
504  }
505  }
506 
507  this->requested = requestedPointcount;
508 
509  unsigned int slice = UINT_MAX/(this->worldSize);
510  unsigned int stateBegin = myRank * slice;
511 
512 //#ifdef HAVE_ZOLTAN2_OMP
513 //#pragma omp parallel for
514 //#endif
515 #ifdef HAVE_ZOLTAN2_OMP
516 #pragma omp parallel
517 #endif
518  {
519  int me = 0;
520  int tsize = 1;
521 #ifdef HAVE_ZOLTAN2_OMP
522  me = omp_get_thread_num();
523  tsize = omp_get_num_threads();
524 #endif
525  unsigned int state = stateBegin + me * slice/(tsize);
526 
527 #ifdef HAVE_ZOLTAN2_OMP
528 #pragma omp for
529 #endif
530  for(lno_t cnt = 0; cnt < requestedPointcount; ++cnt){
531  lno_t iteration = 0;
532  while(1){
533  if(++iteration > MAX_ITER_ALLOWED) {
534  throw "Max number of Iteration is reached for point creation. Check the area criteria or hole coordinates.";
535  }
536  CoordinatePoint <T> p = this->getPoint( this->assignedPrevious + cnt, state);
537 
538  bool isInHole = false;
539  for(lno_t i = 0; i < holeCount; ++i){
540  if(holes[i][0].isInArea(p)){
541  isInHole = true;
542  break;
543  }
544  }
545  if(isInHole) continue;
546  points[cnt].x = p.x;
547 
548  points[cnt].y = p.y;
549  points[cnt].z = p.z;
550  break;
551  }
552  }
553  }
554 //#pragma omp parallel
555  /*
556  {
557 
558  lno_t cnt = 0;
559  lno_t iteration = 0;
560  while (cnt < requestedPointcount){
561  if(++iteration > MAX_ITER_ALLOWED) {
562  throw "Max number of Iteration is reached for point creation. Check the area criteria or hole coordinates.";
563  }
564  CoordinatePoint <T> p = this->getPoint();
565 
566  bool isInHole = false;
567  for(lno_t i = 0; i < holeCount; ++i){
568  if(holes[i][0].isInArea(p)){
569  isInHole = true;
570  break;
571  }
572  }
573  if(isInHole) continue;
574  iteration = 0;
575  points[cnt].x = p.x;
576  points[cnt].y = p.y;
577  points[cnt].z = p.z;
578  ++cnt;
579  }
580  }
581  */
582  }
583 
584  void GetPoints(lno_t requestedPointcount, T **coords/*preallocated sized numPoints*/, lno_t tindex,
585  Hole<T> **holes, lno_t holeCount,
586  float *sharedRatios_, int myRank){
587 
588  for (int i = 0; i < myRank; ++i){
589  //cout << "me:" << myRank << " i:" << i << " s:" << sharedRatios_[i]<< endl;
590  this->assignedPrevious += lno_t(sharedRatios_[i] * this->numPoints);
591  if (gno_t(i) < this->numPoints % this->worldSize){
592  this->assignedPrevious += 1;
593  }
594  }
595 
596  this->requested = requestedPointcount;
597 
598  unsigned int slice = UINT_MAX/(this->worldSize);
599  unsigned int stateBegin = myRank * slice;
600 #ifdef HAVE_ZOLTAN2_OMP
601 #pragma omp parallel
602 #endif
603  {
604  int me = 0;
605  int tsize = 1;
606 #ifdef HAVE_ZOLTAN2_OMP
607  me = omp_get_thread_num();
608  tsize = omp_get_num_threads();
609 #endif
610  unsigned int state = stateBegin + me * (slice/(tsize));
611  /*
612 #pragma omp critical
613  {
614 
615  cout << "myRank:" << me << " stateBeg:" << stateBegin << " tsize:" << tsize << " state:" << state << " slice: " << slice / tsize << endl;
616  }
617  */
618 #ifdef HAVE_ZOLTAN2_OMP
619 #pragma omp for
620 #endif
621  for(lno_t cnt = 0; cnt < requestedPointcount; ++cnt){
622  lno_t iteration = 0;
623  while(1){
624  if(++iteration > MAX_ITER_ALLOWED) {
625  throw "Max number of Iteration is reached for point creation. Check the area criteria or hole coordinates.";
626  }
627  CoordinatePoint <T> p = this->getPoint( this->assignedPrevious + cnt, state);
628 
629  bool isInHole = false;
630  for(lno_t i = 0; i < holeCount; ++i){
631  if(holes[i][0].isInArea(p)){
632  isInHole = true;
633  break;
634  }
635  }
636  if(isInHole) continue;
637  coords[0][cnt + tindex] = p.x;
638  if(this->dimension > 1){
639  coords[1][cnt + tindex] = p.y;
640  if(this->dimension > 2){
641  coords[2][cnt + tindex] = p.z;
642  }
643  }
644  break;
645  }
646  }
647  }
648  }
649 };
650 
651 template <typename T, typename lno_t, typename gno_t>
653 public:
658 
659 
660  virtual T getXCenter(){
661  return center.x;
662  }
663  virtual T getXRadius(){
664  return standartDevx;
665  }
666 
667  CoordinateNormalDistribution(gno_t np_, int dim, CoordinatePoint<T> center_ ,
668  T sd_x, T sd_y, T sd_z, int wSize) :
669  CoordinateDistribution<T,lno_t,gno_t>(np_,dim,wSize),
670  standartDevx(sd_x), standartDevy(sd_y), standartDevz(sd_z)
671  {
672  this->center.x = center_.x;
673  this->center.y = center_.y;
674  this->center.z = center_.z;
675  }
676 
677  virtual CoordinatePoint<T> getPoint(gno_t pindex, unsigned int &state){
678 
679  pindex = 0; // not used in normal distribution.
681 
682  for(int i = 0; i < this->dimension; ++i){
683  switch(i){
684  case 0:
685  p.x = normalDist(this->center.x, this->standartDevx, state);
686  break;
687  case 1:
688  p.y = normalDist(this->center.y, this->standartDevy, state);
689  break;
690  case 2:
691  p.z = normalDist(this->center.z, this->standartDevz, state);
692  break;
693  default:
694  throw "unsupported dimension";
695  }
696  }
697  return p;
698  }
699 
701 private:
702  T normalDist(T center_, T sd, unsigned int &state) {
703  static bool derived=false;
704  static T storedDerivation;
705  T polarsqrt, normalsquared, normal1, normal2;
706  if (!derived) {
707  do {
708  normal1=2.0*( T(rand_r(&state))/T(RAND_MAX) ) - 1.0;
709  normal2=2.0*( T(rand_r(&state))/T(RAND_MAX) ) - 1.0;
710  normalsquared=normal1*normal1+normal2*normal2;
711  } while ( normalsquared>=1.0 || normalsquared == 0.0);
712 
713  polarsqrt=sqrt(-2.0*log(normalsquared)/normalsquared);
714  storedDerivation=normal1*polarsqrt;
715  derived=true;
716  return normal2*polarsqrt*sd + center_;
717  }
718  else {
719  derived=false;
720  return storedDerivation*sd + center_;
721  }
722  }
723 };
724 
725 template <typename T, typename lno_t, typename gno_t>
727 public:
734 
735 
736  virtual T getXCenter(){
737  return (rightMostx - leftMostx)/2 + leftMostx;
738  }
739  virtual T getXRadius(){
740  return (rightMostx - leftMostx)/2;
741  }
742 
743 
744  CoordinateUniformDistribution(gno_t np_, int dim, T l_x, T r_x, T l_y, T r_y,
745  T l_z, T r_z, int wSize ) :
746  CoordinateDistribution<T,lno_t,gno_t>(np_,dim,wSize),
747  leftMostx(l_x), rightMostx(r_x), leftMosty(l_y), rightMosty(r_y),
748  leftMostz(l_z), rightMostz(r_z){}
749 
751  virtual CoordinatePoint<T> getPoint(gno_t pindex, unsigned int &state){
752 
753 
754  pindex = 0; //not used in uniform dist.
756  for(int i = 0; i < this->dimension; ++i){
757  switch(i){
758  case 0:
759  p.x = uniformDist(this->leftMostx, this->rightMostx, state);
760  break;
761  case 1:
762  p.y = uniformDist(this->leftMosty, this->rightMosty, state);
763  break;
764  case 2:
765  p.z = uniformDist(this->leftMostz, this->rightMostz, state);
766  break;
767  default:
768  throw "unsupported dimension";
769  }
770  }
771  return p;
772  }
773 
774 private:
775 
776  T uniformDist(T a, T b, unsigned int &state)
777  {
778  return (b-a)*(rand_r(&state) / double(RAND_MAX)) + a;
779  }
780 };
781 
782 template <typename T, typename lno_t, typename gno_t>
784 public:
791  gno_t along_X, along_Y, along_Z;
792  //T currentX, currentY, currentZ;
794  int myRank;
795  T xstep, ystep, zstep;
796  gno_t xshift, yshift, zshift;
797 
798  virtual T getXCenter(){
799  return (rightMostx - leftMostx)/2 + leftMostx;
800  }
801  virtual T getXRadius(){
802  return (rightMostx - leftMostx)/2;
803  }
804 
805 
806  CoordinateGridDistribution(gno_t alongX, gno_t alongY, gno_t alongZ, int dim,
807  T l_x, T r_x, T l_y, T r_y, T l_z, T r_z ,
808  int myRank_, int wSize) :
809  CoordinateDistribution<T,lno_t,gno_t>(alongX * alongY * alongZ,dim,wSize),
810  leftMostx(l_x), rightMostx(r_x), leftMosty(l_y), rightMosty(r_y), leftMostz(l_z), rightMostz(r_z), myRank(myRank_){
811  //currentX = leftMostx, currentY = leftMosty, currentZ = leftMostz;
812  this->processCnt = 0;
813  this->along_X = alongX; this->along_Y = alongY; this->along_Z = alongZ;
814 
815  if(this->along_X > 1)
816  this->xstep = (rightMostx - leftMostx) / (alongX - 1);
817  else
818  this->xstep = 1;
819  if(this->along_Y > 1)
820  this->ystep = (rightMosty - leftMosty) / (alongY - 1);
821  else
822  this->ystep = 1;
823  if(this->along_Z > 1)
824  this->zstep = (rightMostz - leftMostz) / (alongZ - 1);
825  else
826  this->zstep = 1;
827  xshift = 0; yshift = 0; zshift = 0;
828  }
829 
831  virtual CoordinatePoint<T> getPoint(gno_t pindex, unsigned int &state){
832  //lno_t before = processCnt + this->assignedPrevious;
833  //cout << "before:" << processCnt << " " << this->assignedPrevious << endl;
834  //lno_t xshift = 0, yshift = 0, zshift = 0;
835 
836  //lno_t tmp = before % (this->along_X * this->along_Y);
837  //xshift = tmp % this->along_X;
838  //yshift = tmp / this->along_X;
839  //zshift = before / (this->along_X * this->along_Y);
840 
841  state = 0; //not used here
842  this->zshift = pindex / (along_X * along_Y);
843  this->yshift = (pindex % (along_X * along_Y)) / along_X;
844  this->xshift = (pindex % (along_X * along_Y)) % along_X;
845 
846  //this->xshift = pindex / (along_Z * along_Y);
847  // this->zshift = (pindex % (along_Z * along_Y)) / along_Y;
848  //this->yshift = (pindex % (along_Z * along_Y)) % along_Y;
849 
851  p.x = xshift * this->xstep + leftMostx;
852  p.y = yshift * this->ystep + leftMosty;
853  p.z = zshift * this->zstep + leftMostz;
854 /*
855  ++xshift;
856  if(xshift == this->along_X){
857  ++yshift;
858  xshift = 0;
859  if(yshift == this->along_Y){
860  ++zshift;
861  yshift = 0;
862  }
863  }
864 */
865  /*
866  if(this->processCnt == 0){
867  this->xshift = this->assignedPrevious / (along_Z * along_Y);
868  //this->yshift = (this->assignedPrevious % (along_X * along_Y)) / along_X;
869  this->zshift = (this->assignedPrevious % (along_Z * along_Y)) / along_Y;
870  //this->xshift = (this->assignedPrevious % (along_X * along_Y)) % along_X;
871  this->yshift = (this->assignedPrevious % (along_Z * along_Y)) % along_Y;
872  ++this->processCnt;
873  }
874 
875  CoordinatePoint <T> p;
876  p.x = xshift * this->xstep + leftMostx;
877  p.y = yshift * this->ystep + leftMosty;
878  p.z = zshift * this->zstep + leftMostz;
879 
880  ++yshift;
881  if(yshift == this->along_Y){
882  ++zshift;
883  yshift = 0;
884  if(zshift == this->along_Z){
885  ++xshift;
886  zshift = 0;
887  }
888 
889  }
890  */
891  /*
892  if(this->requested - 1 > this->processCnt){
893  this->processCnt++;
894  } else {
895  cout << "req:" << this->requested << " pr:" <<this->processCnt << endl;
896  cout << "p:" << p.x << " " << p.y << " " << p.z << endl ;
897  cout << "s:" << xshift << " " << yshift << " " << zshift << endl ;
898  cout << "st:" << this->xstep << " " << this->ystep << " " << this->zstep << endl ;
899  }
900  */
901  return p;
902  }
903 
904 private:
905 
906 };
907 
908 template <typename scalar_t, typename lno_t, typename gno_t, typename node_t>
910 private:
911  Hole<scalar_t> **holes; //to represent if there is any hole in the input
912  int holeCount;
913  int coordinate_dimension; //dimension of the geometry
914  gno_t numGlobalCoords; //global number of coordinates requested to be created.
915  lno_t numLocalCoords;
916  float *loadDistributions; //sized as the number of processors, the load of each processor.
917  bool loadDistSet;
918  bool distinctCoordSet;
919  CoordinateDistribution<scalar_t, lno_t,gno_t> **coordinateDistributions;
920  int distributionCount;
921  //CoordinatePoint<scalar_t> *points;
922  scalar_t **coords;
923  scalar_t **wghts;
925  int numWeightsPerCoord;
926  int predistribution;
927  RCP<const Teuchos::Comm<int> > comm;
928  //RCP< Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >tmVector;
929  //RCP< Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >tmwVector;
930  int worldSize;
931  int myRank;
932  scalar_t minx;
933  scalar_t maxx;
934  scalar_t miny;
935  scalar_t maxy;
936  scalar_t minz;
937  scalar_t maxz;
938  std::string outfile;
939  float perturbation_ratio;
940 
941  typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> tMVector_t;
942  typedef Tpetra::Map<lno_t, gno_t, node_t> tMap_t;
943 
944 
945  template <typename tt>
946  tt getParamVal( const Teuchos::ParameterEntry& pe, const std::string &paramname){
947  tt returnVal;
948  try {
949  returnVal = Teuchos::getValue<tt>(pe);
950  }
951  catch (...){
952  throw INVALID(paramname);
953  }
954  return returnVal;
955  }
956 
957  int countChar (std::string inStr, char cntChar){
958  int cnt = 0;
959  for (unsigned int i = 0; i < inStr.size(); ++i){
960  if (inStr[i] == cntChar) {
961  cnt++;
962  }
963  }
964  return cnt;
965  }
966 
967  template <typename tt>
968  std::string toString(tt obj){
969  std::stringstream ss (std::stringstream::in |std::stringstream::out);
970  ss << obj;
971  std::string tmp = "";
972  ss >> tmp;
973  return tmp;
974  }
975 
976  template <typename tt>
977  tt fromString(std::string obj){
978  std::stringstream ss (std::stringstream::in | std::stringstream::out);
979  ss << obj;
980  tt tmp;
981  ss >> tmp;
982  if (ss.fail()){
983  throw "Cannot convert string " + obj;
984  }
985  return tmp;
986  }
987 
988 
989  void splitString(std::string inStr, char splitChar, std::string *outSplittedStr){
990  std::stringstream ss (std::stringstream::in | std::stringstream::out);
991  ss << inStr;
992  int cnt = 0;
993  while (!ss.eof()){
994  std::string tmp = "";
995  std::getline(ss, tmp ,splitChar);
996  outSplittedStr[cnt++] = tmp;
997  }
998  }
999 
1000  /*
1001  void getDistinctCoordinateDescription(std::string distinctDescription){
1002 
1003  this->distinctCoordinates = new bool[this->dimension];
1004  if(distinctDescription != ""){
1005  int argCnt = this->countChar(distinctDescription, ',') + 1;
1006  if(argCnt != this->dimension) {
1007  throw "Invalid parameter count for distinct_coordinates. Size should be equal to dimension.";
1008  }
1009 
1010  std::string *splittedStr = new std::string[argCnt];
1011  splitString(holeDescription, ',', splittedStr);
1012 
1013  for(int i = 0; i < argCnt; ++i){
1014  if(splittedStr[i] == "T"){
1015  distinctCoordinates[i] = true;
1016  } else if(splittedStr[i] == "F"){
1017  distinctCoordinates[i] = false;
1018  } else {
1019  throw "Invalid parameter " + splittedStr[i] + " for distinct_coordinates.";
1020  }
1021  }
1022  delete []splittedStr;
1023  }
1024  else {
1025  for(int i = 0; i < this->dimension; ++i){
1026  distinctCoordinates[i] = false;
1027  }
1028  }
1029  }
1030  */
1031 
1032 
1033  void getCoordinateDistributions(std::string coordinate_distributions){
1034  if(coordinate_distributions == ""){
1035  throw "At least one distribution function must be provided to coordinate generator.";
1036  }
1037 
1038  int argCnt = this->countChar(coordinate_distributions, ',') + 1;
1039  std::string *splittedStr = new std::string[argCnt];
1040  splitString(coordinate_distributions, ',', splittedStr);
1041  coordinateDistributions = (CoordinateDistribution<scalar_t, lno_t,gno_t> **) malloc(sizeof (CoordinateDistribution<scalar_t, lno_t,gno_t> *) * 1);
1042  for(int i = 0; i < argCnt; ){
1043  coordinateDistributions = (CoordinateDistribution<scalar_t, lno_t,gno_t> **)realloc((void *)coordinateDistributions, (this->distributionCount + 1)* sizeof(CoordinateDistribution<scalar_t, lno_t,gno_t> *));
1044 
1045  std::string distName = splittedStr[i++];
1046  gno_t np_ = 0;
1047  if(distName == "NORMAL"){
1048  int reqArg = 5;
1049  if (this->coordinate_dimension == 3){
1050  reqArg = 7;
1051  }
1052  if(i + reqArg > argCnt) {
1053  std::string tmp = toString<int>(reqArg);
1054  throw INVALID_SHAPE_ARG(distName, tmp);
1055  }
1056  np_ = fromString<gno_t>(splittedStr[i++]);
1058 
1059  pp.x = fromString<scalar_t>(splittedStr[i++]);
1060  pp.y = fromString<scalar_t>(splittedStr[i++]);
1061  pp.z = 0;
1062  if(this->coordinate_dimension == 3){
1063  pp.z = fromString<scalar_t>(splittedStr[i++]);
1064  }
1065 
1066  scalar_t sd_x = fromString<scalar_t>(splittedStr[i++]);
1067  scalar_t sd_y = fromString<scalar_t>(splittedStr[i++]);
1068  scalar_t sd_z = 0;
1069  if(this->coordinate_dimension == 3){
1070  sd_z = fromString<scalar_t>(splittedStr[i++]);
1071  }
1072  this->coordinateDistributions[this->distributionCount++] = new CoordinateNormalDistribution<scalar_t, lno_t,gno_t>(np_, this->coordinate_dimension, pp , sd_x, sd_y, sd_z, this->worldSize );
1073 
1074  } else if(distName == "UNIFORM" ){
1075  int reqArg = 5;
1076  if (this->coordinate_dimension == 3){
1077  reqArg = 7;
1078  }
1079  if(i + reqArg > argCnt) {
1080  std::string tmp = toString<int>(reqArg);
1081  throw INVALID_SHAPE_ARG(distName, tmp);
1082  }
1083  np_ = fromString<gno_t>(splittedStr[i++]);
1084  scalar_t l_x = fromString<scalar_t>(splittedStr[i++]);
1085  scalar_t r_x = fromString<scalar_t>(splittedStr[i++]);
1086  scalar_t l_y = fromString<scalar_t>(splittedStr[i++]);
1087  scalar_t r_y = fromString<scalar_t>(splittedStr[i++]);
1088 
1089  scalar_t l_z = 0, r_z = 0;
1090 
1091  if(this->coordinate_dimension == 3){
1092  l_z = fromString<scalar_t>(splittedStr[i++]);
1093  r_z = fromString<scalar_t>(splittedStr[i++]);
1094  }
1095 
1096  this->coordinateDistributions[this->distributionCount++] = new CoordinateUniformDistribution<scalar_t, lno_t,gno_t>( np_, this->coordinate_dimension, l_x, r_x, l_y, r_y, l_z, r_z, this->worldSize );
1097  } else if (distName == "GRID"){
1098  int reqArg = 6;
1099  if(this->coordinate_dimension == 3){
1100  reqArg = 9;
1101  }
1102  if(i + reqArg > argCnt) {
1103  std::string tmp = toString<int>(reqArg);
1104  throw INVALID_SHAPE_ARG(distName, tmp);
1105  }
1106 
1107  gno_t np_x = fromString<gno_t>(splittedStr[i++]);
1108  gno_t np_y = fromString<gno_t>(splittedStr[i++]);
1109  gno_t np_z = 1;
1110 
1111 
1112  if(this->coordinate_dimension == 3){
1113  np_z = fromString<gno_t>(splittedStr[i++]);
1114  }
1115 
1116  np_ = np_x * np_y * np_z;
1117  scalar_t l_x = fromString<scalar_t>(splittedStr[i++]);
1118  scalar_t r_x = fromString<scalar_t>(splittedStr[i++]);
1119  scalar_t l_y = fromString<scalar_t>(splittedStr[i++]);
1120  scalar_t r_y = fromString<scalar_t>(splittedStr[i++]);
1121 
1122  scalar_t l_z = 0, r_z = 0;
1123 
1124  if(this->coordinate_dimension == 3){
1125  l_z = fromString<scalar_t>(splittedStr[i++]);
1126  r_z = fromString<scalar_t>(splittedStr[i++]);
1127  }
1128 
1129  if(np_x < 1 || np_z < 1 || np_y < 1 ){
1130  throw "Provide at least 1 point along each dimension for grid test.";
1131  }
1132  //cout << "ly:" << l_y << " ry:" << r_y << endl;
1133  this->coordinateDistributions[this->distributionCount++] = new CoordinateGridDistribution<scalar_t, lno_t,gno_t>
1134  (np_x, np_y,np_z, this->coordinate_dimension, l_x, r_x,l_y, r_y, l_z, r_z , this->myRank, this->worldSize);
1135 
1136  }
1137  else {
1138  std::string tmp = toString<int>(this->coordinate_dimension);
1139  throw INVALIDSHAPE(distName, tmp);
1140  }
1141  this->numGlobalCoords += (gno_t) np_;
1142  }
1143  delete []splittedStr;
1144  }
1145 
1146  void getProcLoadDistributions(std::string proc_load_distributions){
1147 
1148  this->loadDistributions = new float[this->worldSize];
1149  if(proc_load_distributions == ""){
1150  float r = 1.0 / this->worldSize;
1151  for(int i = 0; i < this->worldSize; ++i){
1152  this->loadDistributions[i] = r;
1153  }
1154  } else{
1155 
1156 
1157  int argCnt = this->countChar(proc_load_distributions, ',') + 1;
1158  if(argCnt != this->worldSize) {
1159  throw "Invalid parameter count load distributions. Given " + toString<int>(argCnt) + " processor size is " + toString<int>(this->worldSize);
1160  }
1161  std::string *splittedStr = new std::string[argCnt];
1162  splitString(proc_load_distributions, ',', splittedStr);
1163  for(int i = 0; i < argCnt; ++i){
1164  this->loadDistributions[i] = fromString<float>(splittedStr[i]);
1165  }
1166  delete []splittedStr;
1167 
1168 
1169  float sum = 0;
1170  for(int i = 0; i < this->worldSize; ++i){
1171  sum += this->loadDistributions[i];
1172  }
1173  if (fabs(sum - 1.0) > 10*std::numeric_limits<float>::epsilon()){
1174  throw "Processor load ratios do not sum to 1.0.";
1175  }
1176  }
1177 
1178  }
1179 
1180  void getHoles(std::string holeDescription){
1181  if(holeDescription == ""){
1182  return;
1183  }
1184  this->holes = (Hole<scalar_t> **) malloc(sizeof (Hole <scalar_t>*));
1185  int argCnt = this->countChar(holeDescription, ',') + 1;
1186  std::string *splittedStr = new std::string[argCnt];
1187  splitString(holeDescription, ',', splittedStr);
1188 
1189  for(int i = 0; i < argCnt; ){
1190  holes = (Hole<scalar_t> **)realloc((void *)holes, (this->holeCount + 1)* sizeof(Hole<scalar_t> *));
1191 
1192  std::string shapeName = splittedStr[i++];
1193  if(shapeName == "SQUARE" && this->coordinate_dimension == 2){
1194  if(i + 3 > argCnt) {
1195  throw INVALID_SHAPE_ARG(shapeName, "3");
1196  }
1198  pp.x = fromString<scalar_t>(splittedStr[i++]);
1199  pp.y = fromString<scalar_t>(splittedStr[i++]);
1200  scalar_t edge = fromString<scalar_t>(splittedStr[i++]);
1201  this->holes[this->holeCount++] = new SquareHole<scalar_t>(pp, edge);
1202  } else if(shapeName == "RECTANGLE" && this->coordinate_dimension == 2){
1203  if(i + 4 > argCnt) {
1204  throw INVALID_SHAPE_ARG(shapeName, "4");
1205  }
1207  pp.x = fromString<scalar_t>(splittedStr[i++]);
1208  pp.y = fromString<scalar_t>(splittedStr[i++]);
1209  scalar_t edgex = fromString<scalar_t>(splittedStr[i++]);
1210  scalar_t edgey = fromString<scalar_t>(splittedStr[i++]);
1211 
1212  this->holes[this->holeCount++] = new RectangleHole<scalar_t>(pp, edgex,edgey);
1213  } else if(shapeName == "CIRCLE" && this->coordinate_dimension == 2){
1214  if(i + 3 > argCnt) {
1215  throw INVALID_SHAPE_ARG(shapeName, "3");
1216  }
1218  pp.x = fromString<scalar_t>(splittedStr[i++]);
1219  pp.y = fromString<scalar_t>(splittedStr[i++]);
1220  scalar_t r = fromString<scalar_t>(splittedStr[i++]);
1221  this->holes[this->holeCount++] = new CircleHole<scalar_t>(pp, r);
1222  } else if(shapeName == "CUBE" && this->coordinate_dimension == 3){
1223  if(i + 4 > argCnt) {
1224  throw INVALID_SHAPE_ARG(shapeName, "4");
1225  }
1227  pp.x = fromString<scalar_t>(splittedStr[i++]);
1228  pp.y = fromString<scalar_t>(splittedStr[i++]);
1229  pp.z = fromString<scalar_t>(splittedStr[i++]);
1230  scalar_t edge = fromString<scalar_t>(splittedStr[i++]);
1231  this->holes[this->holeCount++] = new CubeHole<scalar_t>(pp, edge);
1232  } else if(shapeName == "RECTANGULAR_PRISM" && this->coordinate_dimension == 3){
1233  if(i + 6 > argCnt) {
1234  throw INVALID_SHAPE_ARG(shapeName, "6");
1235  }
1237  pp.x = fromString<scalar_t>(splittedStr[i++]);
1238  pp.y = fromString<scalar_t>(splittedStr[i++]);
1239  pp.z = fromString<scalar_t>(splittedStr[i++]);
1240  scalar_t edgex = fromString<scalar_t>(splittedStr[i++]);
1241  scalar_t edgey = fromString<scalar_t>(splittedStr[i++]);
1242  scalar_t edgez = fromString<scalar_t>(splittedStr[i++]);
1243  this->holes[this->holeCount++] = new RectangularPrismHole<scalar_t>(pp, edgex, edgey, edgez);
1244 
1245  } else if(shapeName == "SPHERE" && this->coordinate_dimension == 3){
1246  if(i + 4 > argCnt) {
1247  throw INVALID_SHAPE_ARG(shapeName, "4");
1248  }
1250  pp.x = fromString<scalar_t>(splittedStr[i++]);
1251  pp.y = fromString<scalar_t>(splittedStr[i++]);
1252  pp.z = fromString<scalar_t>(splittedStr[i++]);
1253  scalar_t r = fromString<scalar_t>(splittedStr[i++]);
1254  this->holes[this->holeCount++] = new SphereHole<scalar_t>(pp, r);
1255  } else {
1256  std::string tmp = toString<int>(this->coordinate_dimension);
1257  throw INVALIDSHAPE(shapeName, tmp);
1258  }
1259  }
1260  delete [] splittedStr;
1261  }
1262 
1263  void getWeightDistribution(std::string *weight_distribution_arr, int wdimension){
1264  int wcount = 0;
1265 
1266  this->wd = new WeightDistribution<scalar_t,scalar_t> *[wdimension];
1267  for(int ii = 0; ii < MAX_WEIGHT_DIM; ++ii){
1268  std::string weight_distribution = weight_distribution_arr[ii];
1269  if(weight_distribution == "") continue;
1270  if(wcount == wdimension) {
1271  throw "Weight Dimension is provided as " + toString<int>(wdimension) + ". More weight distribution is provided.";
1272  }
1273 
1274  int count = this->countChar(weight_distribution, ' ');
1275  std::string *splittedStr = new string[count + 1];
1276  this->splitString(weight_distribution, ' ', splittedStr);
1277  //cout << count << endl;
1278  scalar_t c=1;
1279  scalar_t a1=0,a2=0,a3=0;
1280  scalar_t x1=0,y1=0,z1=0;
1281  scalar_t b1=1,b2=1,b3=1;
1282  scalar_t *steps = NULL;
1283  scalar_t *values= NULL;
1284  int stepCount = 0;
1285  int valueCount = 1;
1286 
1287  for (int i = 1; i < count + 1; ++i){
1288  int assignmentCount = this->countChar(splittedStr[i], '=');
1289  if (assignmentCount != 1){
1290  throw "Error at the argument " + splittedStr[i];
1291  }
1292  std::string parameter_vs_val[2];
1293  this->splitString(splittedStr[i], '=', parameter_vs_val);
1294  std::string parameter = parameter_vs_val[0];
1295  std::string value = parameter_vs_val[1];
1296  //cout << "parameter:" << parameter << " value:" << value << endl;
1297 
1298  if (parameter == "a1"){
1299  a1 = this->fromString<scalar_t>(value);
1300  }
1301  else if (parameter == "a2"){
1302  if(this->coordinate_dimension > 1){
1303  a2 = this->fromString<scalar_t>(value);
1304  }
1305  else {
1306  throw parameter+ " argument is not valid when dimension is " + toString<int>(this->coordinate_dimension);
1307  }
1308 
1309  }
1310  else if (parameter == "a3"){
1311  if(this->coordinate_dimension > 2){
1312  a3 = this->fromString<scalar_t>(value);
1313  }
1314  else {
1315  throw parameter+ " argument is not valid when dimension is " + toString<int>(this->coordinate_dimension);
1316  }
1317  }
1318  else if (parameter == "b1"){
1319  b1 = this->fromString<scalar_t>(value);
1320  }
1321  else if (parameter == "b2"){
1322  if(this->coordinate_dimension > 1){
1323  b2 = this->fromString<scalar_t>(value);
1324  }
1325  else {
1326  throw parameter+ " argument is not valid when dimension is " + toString<int>(this->coordinate_dimension);
1327  }
1328  }
1329  else if (parameter == "b3"){
1330 
1331  if(this->coordinate_dimension > 2){
1332  b3 = this->fromString<scalar_t>(value);
1333  }
1334  else {
1335  throw parameter+ " argument is not valid when dimension is " + toString<int>(this->coordinate_dimension);
1336  }
1337  }
1338  else if (parameter == "c"){
1339  c = this->fromString<scalar_t>(value);
1340  }
1341  else if (parameter == "x1"){
1342  x1 = this->fromString<scalar_t>(value);
1343  }
1344  else if (parameter == "y1"){
1345  if(this->coordinate_dimension > 1){
1346  y1 = this->fromString<scalar_t>(value);
1347  }
1348  else {
1349  throw parameter+ " argument is not valid when dimension is " + toString<int>(this->coordinate_dimension);
1350  }
1351  }
1352  else if (parameter == "z1"){
1353  if(this->coordinate_dimension > 2){
1354  z1 = this->fromString<scalar_t>(value);
1355  }
1356  else {
1357  throw parameter+ " argument is not valid when dimension is " + toString<int>(this->coordinate_dimension);
1358  }
1359  }
1360  else if (parameter == "steps"){
1361  stepCount = this->countChar(value, ',') + 1;
1362  std::string *stepstr = new std::string[stepCount];
1363  this->splitString(value, ',', stepstr);
1364  steps = new scalar_t[stepCount];
1365  for (int j = 0; j < stepCount; ++j){
1366  steps[j] = this->fromString<scalar_t>(stepstr[j]);
1367  }
1368  delete [] stepstr;
1369  }
1370  else if (parameter == "values"){
1371  valueCount = this->countChar(value, ',') + 1;
1372  std::string *stepstr = new std::string[valueCount];
1373  this->splitString(value, ',', stepstr);
1374  values = new scalar_t[valueCount];
1375  for (int j = 0; j < valueCount; ++j){
1376  values[j] = this->fromString<scalar_t>(stepstr[j]);
1377  }
1378  delete [] stepstr;
1379  }
1380  else {
1381  throw "Invalid parameter name at " + splittedStr[i];
1382  }
1383  }
1384 
1385  delete []splittedStr;
1386  if(stepCount + 1!= valueCount){
1387  throw "Step count: " + this->toString<int>(stepCount) + " must be 1 less than value count: " + this->toString<int>(valueCount);
1388  }
1389 
1390 
1391  this->wd[wcount++] = new SteppedEquation<scalar_t,scalar_t>(a1, a2, a3, b1, b2, b3, c, x1, y1, z1, steps, values, stepCount);
1392 
1393  if(stepCount > 0){
1394  delete [] steps;
1395  delete [] values;
1396 
1397  }
1398  }
1399  if(wcount != this->numWeightsPerCoord){
1400  throw "Weight Dimension is provided as " + toString<int>(wdimension) + ". But " + toString<int>(wcount)+" weight distributions are provided.";
1401  }
1402  }
1403 
1404  void parseParams(const Teuchos::ParameterList & params){
1405  try {
1406  std::string holeDescription = "";
1407  std::string proc_load_distributions = "";
1408  std::string distinctDescription = "";
1409  std::string coordinate_distributions = "";
1410  std::string numWeightsPerCoord_parameters[MAX_WEIGHT_DIM];
1411  for (int i = 0; i < MAX_WEIGHT_DIM; ++i){
1412  numWeightsPerCoord_parameters[i] = "";
1413  }
1414 
1415 
1416  for (Teuchos::ParameterList::ConstIterator pit = params.begin(); pit != params.end(); ++pit ){
1417  const std::string &paramName = params.name(pit);
1418 
1419  const Teuchos::ParameterEntry &pe = params.entry(pit);
1420 
1421  if(paramName.find("hole-") == 0){
1422  if(holeDescription == ""){
1423  holeDescription = getParamVal<std::string>(pe, paramName);
1424  }
1425  else {
1426  holeDescription +=","+getParamVal<std::string>(pe, paramName);
1427  }
1428  }
1429  else if(paramName.find("distribution-") == 0){
1430  if(coordinate_distributions == "")
1431  coordinate_distributions = getParamVal<std::string>(pe, paramName);
1432  else
1433  coordinate_distributions += ","+getParamVal<std::string>(pe, paramName);
1434  //cout << coordinate_distributions << endl;
1435  //TODO coordinate distribution description
1436  }
1437 
1438  else if (paramName.find(weight_distribution_string) == 0){
1439  std::string weight_dist_param = paramName.substr(weight_distribution_string.size());
1440  int dash_pos = weight_dist_param.find("-");
1441  std::string distribution_index_string = weight_dist_param.substr(0, dash_pos);
1442  int distribution_index = fromString<int>(distribution_index_string);
1443 
1444  if(distribution_index >= MAX_WEIGHT_DIM){
1445  throw "Given distribution index:" + distribution_index_string + " larger than maximum allowed number of weights:" + toString<int>(MAX_WEIGHT_DIM);
1446  }
1447  numWeightsPerCoord_parameters[distribution_index] += " " + weight_dist_param.substr(dash_pos + 1)+ "="+ getParamVal<std::string>(pe, paramName);
1448  }
1449  else if(paramName == "dim"){
1450  int dim = fromString<int>(getParamVal<std::string>(pe, paramName));
1451  if(dim < 2 && dim > 3){
1452  throw INVALID(paramName);
1453  } else {
1454  this->coordinate_dimension = dim;
1455  }
1456  }
1457  else if(paramName == "wdim"){
1458  int dim = fromString<int>(getParamVal<std::string>(pe, paramName));
1459  if(dim < 1 && dim > MAX_WEIGHT_DIM){
1460  throw INVALID(paramName);
1461  } else {
1462  this->numWeightsPerCoord = dim;
1463  }
1464  }
1465  else if(paramName == "predistribution"){
1466  int pre_distribution = fromString<int>(getParamVal<std::string>(pe, paramName));
1467  if(pre_distribution < 0 && pre_distribution > 3){
1468  throw INVALID(paramName);
1469  } else {
1470  this->predistribution = pre_distribution;
1471  }
1472  }
1473  else if(paramName == "perturbation_ratio"){
1474  float perturbation = fromString<float>(getParamVal<std::string>(pe, paramName));
1475  if(perturbation < 0 && perturbation > 1){
1476  throw INVALID(paramName);
1477  } else {
1478  this->perturbation_ratio = perturbation;
1479  }
1480  }
1481 
1482 
1483  else if(paramName == "proc_load_distributions"){
1484  proc_load_distributions = getParamVal<std::string>(pe, paramName);
1485  this->loadDistSet = true;
1486  }
1487 
1488  else if(paramName == "distinct_coordinates"){
1489  distinctDescription = getParamVal<std::string>(pe, paramName);
1490  if(distinctDescription == "T"){
1491  this->distinctCoordSet = true;
1492  } else if(distinctDescription == "F"){
1493  this->distinctCoordSet = true;
1494  } else {
1495  throw "Invalid parameter for distinct_coordinates: " + distinctDescription + ". Candidates: T or F." ;
1496  }
1497  }
1498 
1499  else if(paramName == "out_file"){
1500  this->outfile = getParamVal<std::string>(pe, paramName);
1501  }
1502 
1503  else {
1504  throw INVALID(paramName);
1505  }
1506  }
1507 
1508 
1509  if(this->coordinate_dimension == 0){
1510  throw "Dimension must be provided to coordinate generator.";
1511  }
1512  /*
1513  if(this->globalNumCoords == 0){
1514  throw "Number of coordinates must be provided to coordinate generator.";
1515  }
1516  */
1517  /*
1518  if(maxx <= minx ){
1519  throw "Error: maxx= "+ toString<scalar_t>(maxx)+ " and minx=" + toString<scalar_t>(minx);
1520  }
1521  if(maxy <= miny ){
1522  throw "Error: maxy= "+ toString<scalar_t>(maxy)+ " and miny=" + toString<scalar_t>(miny);
1523 
1524  }
1525  if(this->dimension == 3 && maxz <= minz ){
1526  throw "Error: maxz= "+ toString<scalar_t>(maxz)+ " and minz=" + toString<scalar_t>(minz);
1527  }
1528  */
1529  if (this->loadDistSet && this->distinctCoordSet){
1530  throw "distinct_coordinates and proc_load_distributions parameters cannot be satisfied together.";
1531  }
1532  this->getHoles(holeDescription);
1533  //this->getDistinctCoordinateDescription(distinctDescription);
1534  this->getProcLoadDistributions(proc_load_distributions);
1535  this->getCoordinateDistributions(coordinate_distributions);
1536  this->getWeightDistribution(numWeightsPerCoord_parameters, this->numWeightsPerCoord);
1537  /*
1538  if(this->numGlobalCoords <= 0){
1539  throw "Must have at least 1 point";
1540  }
1541  */
1542  }
1543  catch(std::string s){
1544  throw s;
1545  }
1546 
1547  catch(char * s){
1548  throw s;
1549  }
1550 
1551  catch(char const* s){
1552  throw s;
1553  }
1554 
1555  }
1556 public:
1557 
1559  if(this->holes){
1560  for (int i = 0; i < this->holeCount; ++i){
1561  delete this->holes[i];
1562  }
1563  free (this->holes);
1564  }
1565 
1566 
1567  delete []loadDistributions; //sized as the number of processors, the load of each processor.
1568  //delete []distinctCoordinates; // if processors have different or same range for coordinates to be created.
1569  if(coordinateDistributions){
1570 
1571  for (int i = 0; i < this->distributionCount; ++i){
1572  delete this->coordinateDistributions[i];
1573  }
1574  free (this->coordinateDistributions);
1575  }
1576  if (this->wd){
1577  for (int i = 0; i < this->numWeightsPerCoord; ++i){
1578  delete this->wd[i];
1579  }
1580  delete []this->wd;
1581  }
1582 
1583  if(this->numWeightsPerCoord){
1584  for(int i = 0; i < this->numWeightsPerCoord; ++i)
1585  delete [] this->wghts[i];
1586  delete []this->wghts;
1587  }
1588  if(this->coordinate_dimension){
1589  for(int i = 0; i < this->coordinate_dimension; ++i)
1590  delete [] this->coords[i];
1591  delete [] this->coords;
1592  }
1593  //delete []this->points;
1594  }
1595 
1597  cout <<"\nGeometric Generator Parameter File Format:" << endl;
1598  cout <<"- dim=Coordinate Dimension: 2 or 3" << endl;
1599  cout <<"- Available distributions:" << endl;
1600  cout <<"\tUNIFORM: -> distribution-1=UNIFORM,NUMCOORDINATES,XMIN,XMAX,YMIN,YMAX{,ZMIN,ZMAX}" << endl;
1601  cout <<"\tGRID: -> distribution-2=GRID,XLENGTH,YLENGTH{,ZLENGTH},XMIN,XMAX,YMIN,YMAX{,ZMIN,ZMAX}" << endl;
1602  cout <<"\tNORMAL: -> distribution-3=NORMAL,XCENTER,YCENTER{,ZCENTER},XSD,YSD,{,ZSD}" << endl;
1603  cout <<"- wdim=numWeightsPerCoord: There should be as many weight function as number of weights per coord." << endl;
1604  cout <<"- Weight Equation: w = (a1 * (x - x1)^b1) + (a2 * (y - y1)^b2) + (a3 * (z - z1)^b3) + c" << endl;
1605  cout << "Parameter settings:" << endl;
1606  cout << "\tWeightDistribution-1-a1=a1 " << endl;
1607  cout << "\tWeightDistribution-1-a2=a2 " << endl;
1608  cout << "\tWeightDistribution-1-a3=a3 " << endl;
1609  cout << "\tWeightDistribution-1-b1=b1 " << endl;
1610  cout << "\tWeightDistribution-1-b2=b2 " << endl;
1611  cout << "\tWeightDistribution-1-b3=b3 " << endl;
1612  cout << "\tWeightDistribution-1-x1=x1 " << endl;
1613  cout << "\tWeightDistribution-1-y1=y1 " << endl;
1614  cout << "\tWeightDistribution-1-z1=z1 " << endl;
1615  cout << "\tWeightDistribution-1-c=c " << endl;
1616  cout << "\tIt is possible to set step function to the result of weight equation." << endl;
1617  cout << "\tWeightDistribution-1-steps=step1,step2,step3:increasing order" << endl;
1618  cout << "\tWeightDistribution-1-values=val1,val2,val3,val4." << endl;
1619  cout << "\t\tIf w < step1 -> w = val1" << endl;
1620  cout << "\t\tElse if w < step2 -> w = val2" << endl;
1621  cout << "\t\tElse if w < step3 -> w = val3" << endl;
1622  cout << "\t\tElse -> w = val4" << endl;
1623  cout <<"- Holes:" << endl;
1624  cout << "\thole-1:SPHERE,XCENTER,YCENTER,ZCENTER,RADIUS (only for dim=3)" << endl;
1625  cout << "\thole-2:CUBE,XCENTER,YCENTER,ZCENTER,EDGE (only for dim=3)" << endl;
1626  cout << "\thole-3:RECTANGULAR_PRISM,XCENTER,YCENTER,ZCENTER,XEDGE,YEDGE,ZEDGE (only for dim=3)" << endl;
1627  cout << "\thole-4:SQUARE,XCENTER,YCENTER,EDGE (only for dim=2)" << endl;
1628  cout << "\thole-5:RECTANGLE,XCENTER,YCENTER,XEDGE,YEDGE (only for dim=2)" << endl;
1629  cout << "\thole-6:CIRCLE,XCENTER,YCENTER,RADIUS (only for dim=2)" << endl;
1630  cout << "- out_file:out_file_path : if provided output will be written to files." << endl;
1631  cout << "- proc_load_distributions:ratio_0, ratio_1, ratio_2....ratio_n. Loads of each processor, should be as many as MPI ranks and should sum up to 1." << endl;
1632  cout << "- predistribution:distribution_option. Predistribution of the coordinates to the processors. 0 for NONE, 1 RCB, 2 MJ, 3 BLOCK." << endl;
1633  cout << "- perturbation_ratio:the percent of the local data which will be perturbed in order to simulate the changes in the dynamic partitioning. Float value between 0 and 1." << endl;
1634  }
1635 
1636  GeometricGenerator(Teuchos::ParameterList &params, const RCP<const Teuchos::Comm<int> > & comm_){
1637  this->wd = NULL;
1638  this->comm = comm_;
1639  this->holes = NULL; //to represent if there is any hole in the input
1640  this->coordinate_dimension = 0; //dimension of the geometry
1641  this->numWeightsPerCoord = 0;
1642  this->worldSize = comm_->getSize(); //comminication world object.
1643  this->numGlobalCoords = 0; //global number of coordinates requested to be created.
1644  this->loadDistributions = NULL; //sized as the number of processors, the load of each processor.
1645  //this->distinctCoordinates = NULL; // if processors have different or same range for coordinates to be created.
1646  this->coordinateDistributions = NULL;
1647  this->holeCount = 0;
1648  this->distributionCount = 0;
1649  this->outfile = "";
1650  this->predistribution = 0;
1651  this->perturbation_ratio = 0;
1652  //this->points = NULL;
1653 
1654  /*
1655  this->minx = 0;
1656  this->maxx = 0;
1657  this->miny = 0;
1658  this->maxy = 0;
1659  this->minz = 0;
1660  this->maxz = 0;
1661  */
1662  this->loadDistSet = false;
1663  this->distinctCoordSet = false;
1664  this->myRank = comm_->getRank();
1665 
1666  try {
1667  this->parseParams(params);
1668  }
1669  catch(std::string s){
1670  if(myRank == 0){
1671  print_description();
1672  }
1673  throw s;
1674  }
1675 
1676 
1677 
1678 
1679  lno_t myPointCount = 0;
1680  this->numGlobalCoords = 0;
1681 
1682  gno_t prefixSum = 0;
1683  for(int i = 0; i < this->distributionCount; ++i){
1684  for(int ii = 0; ii < this->worldSize; ++ii){
1685  lno_t increment = lno_t (this->coordinateDistributions[i]->numPoints * this->loadDistributions[ii]);
1686  if (gno_t(ii) < this->coordinateDistributions[i]->numPoints % this->worldSize){
1687  increment += 1;
1688  }
1689  this->numGlobalCoords += increment;
1690  if(ii < myRank){
1691  prefixSum += increment;
1692  }
1693  }
1694  myPointCount += lno_t(this->coordinateDistributions[i]->numPoints * this->loadDistributions[myRank]);
1695  if (gno_t(myRank) < this->coordinateDistributions[i]->numPoints % this->worldSize){
1696  myPointCount += 1;
1697  }
1698  }
1699 
1700  this->coords = new scalar_t *[this->coordinate_dimension];
1701  for(int i = 0; i < this->coordinate_dimension; ++i){
1702  this->coords[i] = new scalar_t[myPointCount];
1703  }
1704 
1705  for (int ii = 0; ii < this->coordinate_dimension; ++ii){
1706 #ifdef HAVE_ZOLTAN2_OMP
1707 #pragma omp parallel for
1708 #endif
1709  for(lno_t i = 0; i < myPointCount; ++i){
1710  this->coords[ii][i] = 0;
1711  }
1712  }
1713 
1714  this->numLocalCoords = 0;
1715  srand ((myRank + 1) * this->numLocalCoords);
1716  for (int i = 0; i < distributionCount; ++i){
1717 
1718  lno_t requestedPointCount = lno_t(this->coordinateDistributions[i]->numPoints * this->loadDistributions[myRank]);
1719  if (gno_t(myRank) < this->coordinateDistributions[i]->numPoints % this->worldSize){
1720  requestedPointCount += 1;
1721  }
1722  //cout << "req:" << requestedPointCount << endl;
1723  //this->coordinateDistributions[i]->GetPoints(requestedPointCount,this->points + this->numLocalCoords, this->holes, this->holeCount, this->loadDistributions, myRank);
1724  this->coordinateDistributions[i]->GetPoints(requestedPointCount,this->coords, this->numLocalCoords, this->holes, this->holeCount, this->loadDistributions, myRank);
1725  this->numLocalCoords += requestedPointCount;
1726  }
1727 
1728  /*
1729 
1730  if (this->myRank >= 0){
1731  for(lno_t i = 0; i < this->numLocalCoords; ++i){
1732 
1733  cout <<"me:" << this->myRank << " "<< this->coords[0][i];
1734  if(this->coordinate_dimension > 1){
1735  cout << " " << this->coords[1][i];
1736  }
1737  if(this->coordinate_dimension > 2){
1738  cout << " " << this->coords[2][i];
1739  }
1740  cout << std::endl;
1741  }
1742  }
1743  */
1744 
1745 
1746 
1747  if (this->predistribution){
1748  redistribute();
1749  }
1750 
1751 
1752 
1753  int scale = 3;
1754  if (this->perturbation_ratio > 0){
1755  this->perturb_data(this->perturbation_ratio, scale);
1756  }
1757  /*
1758  if (this->myRank >= 0){
1759  cout << "after distribution" << endl;
1760  for(lno_t i = 0; i < this->numLocalCoords; ++i){
1761 
1762  cout <<"me:" << this->myRank << " " << this->coords[0][i];
1763  if(this->coordinate_dimension > 1){
1764  cout << " " << this->coords[1][i];
1765  }
1766  if(this->coordinate_dimension > 2){
1767  cout << " " << this->coords[2][i];
1768  }
1769  cout << std::endl;
1770  }
1771  }
1772 
1773  */
1774 
1775 
1776  if (this->distinctCoordSet){
1777  //TODO: Partition and migration.
1778  }
1779 
1780  if(this->outfile != ""){
1781 
1782  std::ofstream myfile;
1783  myfile.open ((this->outfile + toString<int>(myRank)).c_str());
1784  for(lno_t i = 0; i < this->numLocalCoords; ++i){
1785 
1786  myfile << this->coords[0][i];
1787  if(this->coordinate_dimension > 1){
1788  myfile << " " << this->coords[1][i];
1789  }
1790  if(this->coordinate_dimension > 2){
1791  myfile << " " << this->coords[2][i];
1792  }
1793  myfile << std::endl;
1794  }
1795  myfile.close();
1796 
1797  if (this->myRank == 0){
1798  std::ofstream gnuplotfile("gnu.gnuplot");
1799  for(int i = 0; i < this->worldSize; ++i){
1800  string s = "splot";
1801  if (this->coordinate_dimension == 2){
1802  s = "plot";
1803  }
1804  if (i > 0){
1805  s = "replot";
1806  }
1807  gnuplotfile << s << " \"" << (this->outfile + toString<int>(i)) << "\"" << endl;
1808  }
1809  gnuplotfile << "pause -1" << endl;
1810  gnuplotfile.close();
1811  }
1812  }
1813 
1814 
1815 
1816  /*
1817  Zoltan2::XpetraMultiVectorAdapter < Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> > xmv (RCP < Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> > (tmVector));
1818 
1819  RCP< Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >tmVector2;
1820  Zoltan2::PartitioningSolution< Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> > solution;
1821  xmv.applyPartitioningSolution<Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >(this->tmVector, &tmVector2, solution);
1822  */
1823  if (this->numWeightsPerCoord > 0){
1824  this->wghts = new scalar_t *[this->numWeightsPerCoord];
1825  for(int i = 0; i < this->numWeightsPerCoord; ++i){
1826  this->wghts[i] = new scalar_t[this->numLocalCoords];
1827  }
1828  }
1829 
1830  for(int ii = 0; ii < this->numWeightsPerCoord; ++ii){
1831  switch(this->coordinate_dimension){
1832  case 1:
1833 #ifdef HAVE_ZOLTAN2_OMP
1834 #pragma omp parallel for
1835 #endif
1836  for (lno_t i = 0; i < this->numLocalCoords; ++i){
1837  this->wghts[ii][i] = this->wd[ii]->get1DWeight(this->coords[0][i]);
1838  }
1839  break;
1840  case 2:
1841 #ifdef HAVE_ZOLTAN2_OMP
1842 #pragma omp parallel for
1843 #endif
1844  for (lno_t i = 0; i < this->numLocalCoords; ++i){
1845  this->wghts[ii][i] = this->wd[ii]->get2DWeight(this->coords[0][i], this->coords[1][i]);
1846  }
1847  break;
1848  case 3:
1849 #ifdef HAVE_ZOLTAN2_OMP
1850 #pragma omp parallel for
1851 #endif
1852  for (lno_t i = 0; i < this->numLocalCoords; ++i){
1853  this->wghts[ii][i] = this->wd[ii]->get3DWeight(this->coords[0][i], this->coords[1][i], this->coords[2][i]);
1854  }
1855  break;
1856  }
1857  }
1858  }
1859 
1860  //############################################################//
1862  //############################################################//
1863  void perturb_data(float used_perturbation_ratio, int scale){
1864  scalar_t *maxCoords= new scalar_t [this->coordinate_dimension];
1865  scalar_t *minCoords= new scalar_t [this->coordinate_dimension];
1866  for (int dim = 0; dim < this->coordinate_dimension; ++dim){
1867  minCoords[dim] = maxCoords[dim] = this->coords[dim][0];
1868  for (lno_t i = 1; i < this->numLocalCoords; ++i){
1869  if (minCoords[dim] > this->coords[dim][i]){
1870  minCoords[dim] = this->coords[dim][i];
1871  }
1872 
1873  if (maxCoords[dim] < this->coords[dim][i]){
1874  maxCoords[dim] = this->coords[dim][i];
1875  }
1876  }
1877 
1878 
1879 
1880 
1881  scalar_t center = (maxCoords[dim] + minCoords[dim]) / 2;
1882 
1883  minCoords[dim] = center - (center - minCoords[dim]) * scale;
1884  maxCoords[dim] = (maxCoords[dim] - center) * scale + center;
1885 
1886  }
1887 
1888  gno_t numLocalToPerturb = gno_t(this->numLocalCoords*used_perturbation_ratio);
1889  //cout << "numLocalToPerturb :" << numLocalToPerturb << endl;
1890  for (int dim = 0; dim < this->coordinate_dimension; ++dim){
1891  scalar_t range = maxCoords[dim] - minCoords[dim];
1892  for (gno_t i = 0; i < numLocalToPerturb; ++i){
1893  this->coords[dim][i] = (rand() / double (RAND_MAX)) * (range) + minCoords[dim];
1894 
1895  }
1896  }
1897  delete []maxCoords;
1898  delete []minCoords;
1899  }
1900 
1901  //############################################################//
1903  //############################################################//
1904 
1905  //Returns the partitioning dimension as even as possible
1906  void getBestSurface (int remaining, int *dimProcs, int dim, int currentDim, int &bestSurface, int *bestDimProcs){
1907 
1908  if (currentDim < dim - 1){
1909  int ipx = 1;
1910  while(ipx <= remaining) {
1911  if(remaining % ipx == 0) {
1912  int nremain = remaining / ipx;
1913  dimProcs[currentDim] = ipx;
1914  getBestSurface (nremain, dimProcs, dim, currentDim + 1, bestSurface, bestDimProcs);
1915  }
1916  ipx++;
1917  }
1918  }
1919  else {
1920  dimProcs[currentDim] = remaining;
1921  int surface = 0;
1922  for (int i = 0; i < dim; ++i) surface += dimProcs[i];
1923  if (surface < bestSurface){
1924  bestSurface = surface;
1925  for (int i = 0; i < dim; ++i) bestDimProcs[i] = dimProcs[i];
1926  }
1927  }
1928 
1929  }
1930 
1931  //returns min and max coordinates along each dimension
1932  void getMinMaxCoords(scalar_t *globalMaxCoords, scalar_t *globalMinCoords){
1933  scalar_t *maxCoords= new scalar_t [this->coordinate_dimension];
1934  scalar_t *minCoords= new scalar_t [this->coordinate_dimension];
1935  for (int dim = 0; dim < this->coordinate_dimension; ++dim){
1936  minCoords[dim] = maxCoords[dim] = this->coords[dim][0];
1937  for (lno_t i = 1; i < this->numLocalCoords; ++i){
1938  if (minCoords[dim] > this->coords[dim][i]){
1939  minCoords[dim] = this->coords[dim][i];
1940  }
1941 
1942  if (maxCoords[dim] < this->coords[dim][i]){
1943  maxCoords[dim] = this->coords[dim][i];
1944  }
1945  }
1946  }
1947 
1948  reduceAll<int, scalar_t>( *(this->comm), Teuchos::REDUCE_MAX,
1949  this->coordinate_dimension,
1950  maxCoords,
1951  globalMaxCoords);
1952 
1953 
1954  reduceAll<int, scalar_t>( *(this->comm), Teuchos::REDUCE_MIN,
1955  this->coordinate_dimension,
1956  minCoords,
1957  globalMinCoords);
1958 
1959  delete []maxCoords;
1960  delete []minCoords;
1961  }
1962 
1963 
1964  //performs a block partitioning.
1965  //then distributes the points of the overloaded parts to underloaded parts.
1966  void blockPartition(int *coordinate_grid_parts){
1967 
1968 #ifdef _MSC_VER
1969  typedef SSIZE_T ssize_t;
1970 #endif
1971 
1972  //############################################################//
1974  //############################################################//
1975 
1976  scalar_t *maxCoords= new scalar_t [this->coordinate_dimension];
1977  scalar_t *minCoords= new scalar_t [this->coordinate_dimension];
1978  //global min and max coordinates in each dimension.
1979  this->getMinMaxCoords(maxCoords, minCoords);
1980 
1981 
1982  //############################################################//
1984  //############################################################//
1985  int remaining = this->worldSize;
1986  int coord_dim = this->coordinate_dimension;
1987  int *dimProcs = new int[coord_dim];
1988  int *bestDimProcs = new int[coord_dim];
1989  int currentDim = 0;
1990 
1991  int bestSurface = 0;
1992  dimProcs[0] = remaining;
1993  for (int i = 1; i < coord_dim; ++i) dimProcs[i] = 1;
1994  for (int i = 0; i < coord_dim; ++i) bestSurface += dimProcs[i];
1995  for (int i = 0; i < coord_dim; ++i) bestDimProcs[i] = dimProcs[i];
1996  //divides the parts into dimensions as even as possible.
1997  getBestSurface ( remaining, dimProcs, coord_dim, currentDim, bestSurface, bestDimProcs);
1998 
1999 
2000  delete []dimProcs;
2001 
2002  //############################################################//
2004  //############################################################//
2005  int *shiftProcCount = new int[coord_dim];
2006  //how the consecutive parts along a dimension
2007  //differs in the index.
2008 
2009  int remainingProc = this->worldSize;
2010  for (int dim = 0; dim < coord_dim; ++dim){
2011  remainingProc = remainingProc / bestDimProcs[dim];
2012  shiftProcCount[dim] = remainingProc;
2013  }
2014 
2015  scalar_t *dim_slices = new scalar_t[coord_dim];
2016  for (int dim = 0; dim < coord_dim; ++dim){
2017  scalar_t dim_range = maxCoords[dim] - minCoords[dim];
2018  scalar_t slice = dim_range / bestDimProcs[dim];
2019  dim_slices[dim] = slice;
2020  }
2021 
2022  //############################################################//
2024  //############################################################//
2025 
2026  gno_t *numPointsInParts = new gno_t[this->worldSize];
2027  gno_t *numGlobalPointsInParts = new gno_t[this->worldSize];
2028  gno_t *numPointsInPartsInclusiveUptoMyIndex = new gno_t[this->worldSize];
2029 
2030  gno_t *partBegins = new gno_t [this->worldSize];
2031  gno_t *partNext = new gno_t[this->numLocalCoords];
2032 
2033 
2034  for (lno_t i = 0; i < this->numLocalCoords; ++i){
2035  partNext[i] = -1;
2036  }
2037  for (int i = 0; i < this->worldSize; ++i) {
2038  partBegins[i] = 1;
2039  }
2040 
2041  for (int i = 0; i < this->worldSize; ++i)
2042  numPointsInParts[i] = 0;
2043 
2044  for (lno_t i = 0; i < this->numLocalCoords; ++i){
2045  int partIndex = 0;
2046  for (int dim = 0; dim < coord_dim; ++dim){
2047  int shift = int ((this->coords[dim][i] - minCoords[dim]) / dim_slices[dim]);
2048  if (shift >= bestDimProcs[dim]){
2049  shift = bestDimProcs[dim] - 1;
2050  }
2051  shift = shift * shiftProcCount[dim];
2052  partIndex += shift;
2053  }
2054  numPointsInParts[partIndex] += 1;
2055  coordinate_grid_parts[i] = partIndex;
2056 
2057  partNext[i] = partBegins[partIndex];
2058  partBegins[partIndex] = i;
2059 
2060  }
2061 
2062  //############################################################//
2064  //############################################################//
2065  reduceAll<int, gno_t>( *(this->comm), Teuchos::REDUCE_SUM,
2066  this->worldSize,
2067  numPointsInParts,
2068  numGlobalPointsInParts);
2069 
2070 
2071  Teuchos::scan<int,gno_t>(
2072  *(this->comm), Teuchos::REDUCE_SUM,
2073  this->worldSize,
2074  numPointsInParts,
2075  numPointsInPartsInclusiveUptoMyIndex
2076  );
2077 
2078 
2079 
2080 
2081  /*
2082  gno_t totalSize = 0;
2083  for (int i = 0; i < this->worldSize; ++i){
2084  totalSize += numPointsInParts[i];
2085  }
2086  if (totalSize != this->numLocalCoords){
2087  cout << "me:" << this->myRank << " ts:" << totalSize << " nl:" << this->numLocalCoords << endl;
2088  }
2089  */
2090 
2091 
2092  //cout << "me:" << this->myRank << " ilk print" << endl;
2093 
2094  gno_t optimal_num = gno_t(this->numGlobalCoords/double(this->worldSize)+0.5);
2095 #ifdef printparts
2096  if (this->myRank == 0){
2097  gno_t totalSize = 0;
2098  for (int i = 0; i < this->worldSize; ++i){
2099  cout << "me:" << this->myRank << " NumPoints in part:" << i << " is: " << numGlobalPointsInParts[i] << endl;
2100  totalSize += numGlobalPointsInParts[i];
2101  }
2102  cout << "Total:" << totalSize << " ng:" << this->numGlobalCoords << endl;
2103  cout << "optimal_num:" << optimal_num << endl;
2104  }
2105 #endif
2106  ssize_t *extraInPart = new ssize_t [this->worldSize];
2107 
2108  ssize_t extraExchange = 0;
2109  for (int i = 0; i < this->worldSize; ++i){
2110  extraInPart[i] = numGlobalPointsInParts[i] - optimal_num;
2111  extraExchange += extraInPart[i];
2112  }
2113  if (extraExchange != 0){
2114  int addition = -1;
2115  if (extraExchange < 0) addition = 1;
2116  for (ssize_t i = 0; i < extraExchange; ++i){
2117  extraInPart[i % this->worldSize] += addition;
2118  }
2119  }
2120 
2121  //############################################################//
2123  //############################################################//
2124 
2125  int overloadedPartCount = 0;
2126  int *overloadedPartIndices = new int [this->worldSize];
2127 
2128 
2129  int underloadedPartCount = 0;
2130  int *underloadedPartIndices = new int [this->worldSize];
2131 
2132  for (int i = 0; i < this->worldSize; ++i){
2133  if(extraInPart[i] > 0){
2134  overloadedPartIndices[overloadedPartCount++] = i;
2135  }
2136  else if(extraInPart[i] < 0){
2137  underloadedPartIndices[underloadedPartCount++] = i;
2138  }
2139  }
2140 
2141  int underloadpartindex = underloadedPartCount - 1;
2142 
2143 
2144  int numPartsISendFrom = 0;
2145  int *mySendFromParts = new int[this->worldSize * 2];
2146  gno_t *mySendFromPartsCounts = new gno_t[this->worldSize * 2];
2147 
2148  int numPartsISendTo = 0;
2149  int *mySendParts = new int[this->worldSize * 2];
2150  gno_t *mySendCountsToParts = new gno_t[this->worldSize * 2];
2151 
2152 
2153  //############################################################//
2158  //############################################################//
2159  for (int i = overloadedPartCount - 1; i >= 0; --i){
2160  //get the overloaded part
2161  //the overload
2162  int overloadPart = overloadedPartIndices[i];
2163  gno_t overload = extraInPart[overloadPart];
2164  gno_t myload = numPointsInParts[overloadPart];
2165 
2166 
2167  //the inclusive load of the processors up to me
2168  gno_t inclusiveLoadUpToMe = numPointsInPartsInclusiveUptoMyIndex[overloadPart];
2169 
2170  //the exclusive load of the processors up to me
2171  gno_t exclusiveLoadUptoMe = inclusiveLoadUpToMe - myload;
2172 
2173 
2174  if (exclusiveLoadUptoMe >= overload){
2175  //this processor does not have to convert anything.
2176  //for this overloaded part.
2177  //set the extra for this processor to zero.
2178  overloadedPartIndices[i] = -1;
2179  extraInPart[overloadPart] = 0;
2180  //then consume underloaded parts.
2181  while (overload > 0){
2182  int nextUnderLoadedPart = underloadedPartIndices[underloadpartindex];
2183  ssize_t underload = extraInPart[nextUnderLoadedPart];
2184  ssize_t left = overload + underload;
2185 
2186  if(left >= 0){
2187  extraInPart[nextUnderLoadedPart] = 0;
2188  underloadedPartIndices[underloadpartindex--] = -1;
2189  }
2190  else {
2191  extraInPart[nextUnderLoadedPart] = left;
2192  }
2193  overload = left;
2194  }
2195  }
2196  else if (exclusiveLoadUptoMe < overload){
2197  //if the previous processors load is not enough
2198  //then this processor should convert some of its elements.
2199  gno_t mySendCount = overload - exclusiveLoadUptoMe;
2200  //how much more needed.
2201  gno_t sendAfterMe = 0;
2202  //if my load is not enough
2203  //then the next processor will continue to convert.
2204  if (mySendCount > myload){
2205  sendAfterMe = mySendCount - myload;
2206  mySendCount = myload;
2207  }
2208 
2209 
2210  //this processors will convert from overloaded part
2211  //as many as mySendCount items.
2212  mySendFromParts[numPartsISendFrom] = overloadPart;
2213  mySendFromPartsCounts[numPartsISendFrom++] = mySendCount;
2214 
2215  //first consume underloaded parts for the previous processors.
2216  while (exclusiveLoadUptoMe > 0){
2217  int nextUnderLoadedPart = underloadedPartIndices[underloadpartindex];
2218  ssize_t underload = extraInPart[nextUnderLoadedPart];
2219  ssize_t left = exclusiveLoadUptoMe + underload;
2220 
2221  if(left >= 0){
2222  extraInPart[nextUnderLoadedPart] = 0;
2223  underloadedPartIndices[underloadpartindex--] = -1;
2224  }
2225  else {
2226  extraInPart[nextUnderLoadedPart] = left;
2227  }
2228  exclusiveLoadUptoMe = left;
2229  }
2230 
2231  //consume underloaded parts for my load.
2232  while (mySendCount > 0){
2233  int nextUnderLoadedPart = underloadedPartIndices[underloadpartindex];
2234  ssize_t underload = extraInPart[nextUnderLoadedPart];
2235  ssize_t left = mySendCount + underload;
2236 
2237  if(left >= 0){
2238  mySendParts[numPartsISendTo] = nextUnderLoadedPart;
2239  mySendCountsToParts[numPartsISendTo++] = -underload;
2240 
2241  extraInPart[nextUnderLoadedPart] = 0;
2242  underloadedPartIndices[underloadpartindex--] = -1;
2243  }
2244  else {
2245  extraInPart[nextUnderLoadedPart] = left;
2246 
2247  mySendParts[numPartsISendTo] = nextUnderLoadedPart;
2248  mySendCountsToParts[numPartsISendTo++] = mySendCount;
2249 
2250  }
2251  mySendCount = left;
2252  }
2253  //consume underloaded parts for the load of the processors after my index.
2254  while (sendAfterMe > 0){
2255  int nextUnderLoadedPart = underloadedPartIndices[underloadpartindex];
2256  ssize_t underload = extraInPart[nextUnderLoadedPart];
2257  ssize_t left = sendAfterMe + underload;
2258 
2259  if(left >= 0){
2260  extraInPart[nextUnderLoadedPart] = 0;
2261  underloadedPartIndices[underloadpartindex--] = -1;
2262  }
2263  else {
2264  extraInPart[nextUnderLoadedPart] = left;
2265  }
2266  sendAfterMe = left;
2267  }
2268  }
2269  }
2270 
2271 
2272  //############################################################//
2274  //############################################################//
2275  for (int i = 0 ; i < numPartsISendFrom; ++i){
2276 
2277  //get the part from which the elements will be converted.
2278  int sendFromPart = mySendFromParts[i];
2279  ssize_t sendCount = mySendFromPartsCounts[i];
2280  while(sendCount > 0){
2281  int partToSendIndex = numPartsISendTo - 1;
2282  int partToSend = mySendParts[partToSendIndex];
2283 
2284  int sendCountToThisPart = mySendCountsToParts[partToSendIndex];
2285 
2286  //determine which part i should convert to
2287  //and how many to this part.
2288  if (sendCountToThisPart <= sendCount){
2289  mySendParts[partToSendIndex] = 0;
2290  mySendCountsToParts[partToSendIndex] = 0;
2291  --numPartsISendTo;
2292  sendCount -= sendCountToThisPart;
2293  }
2294  else {
2295  mySendCountsToParts[partToSendIndex] = sendCountToThisPart - sendCount;
2296  sendCountToThisPart = sendCount;
2297  sendCount = 0;
2298  }
2299 
2300 
2301  gno_t toChange = partBegins[sendFromPart];
2302  gno_t previous_begin = partBegins[partToSend];
2303 
2304  //do the conversion.
2305  for (int k = 0; k < sendCountToThisPart - 1; ++k){
2306  coordinate_grid_parts[toChange] = partToSend;
2307  toChange = partNext[toChange];
2308  }
2309  coordinate_grid_parts[toChange] = partToSend;
2310 
2311  gno_t newBegin = partNext[toChange];
2312  partNext[toChange] = previous_begin;
2313  partBegins[partToSend] = partBegins[sendFromPart];
2314  partBegins[sendFromPart] = newBegin;
2315  }
2316  }
2317 
2318  //if (this->myRank == 0) cout << "4" << endl;
2319 
2320 #ifdef printparts
2321 
2322 
2323  for (int i = 0; i < this->worldSize; ++i) numPointsInParts[i] = 0;
2324 
2325  for (int i = 0; i < this->numLocalCoords; ++i){
2326  numPointsInParts[coordinate_grid_parts[i]] += 1;
2327  }
2328 
2329  reduceAll<int, gno_t>( *(this->comm), Teuchos::REDUCE_SUM,
2330  this->worldSize,
2331  numPointsInParts,
2332  numGlobalPointsInParts);
2333  if (this->myRank == 0){
2334  cout << "reassigning" << endl;
2335  gno_t totalSize = 0;
2336  for (int i = 0; i < this->worldSize; ++i){
2337  cout << "me:" << this->myRank << " NumPoints in part:" << i << " is: " << numGlobalPointsInParts[i] << endl;
2338  totalSize += numGlobalPointsInParts[i];
2339  }
2340  cout << "Total:" << totalSize << " ng:" << this->numGlobalCoords << endl;
2341  }
2342 #endif
2343  delete []mySendCountsToParts;
2344  delete []mySendParts;
2345  delete []mySendFromPartsCounts;
2346  delete []mySendFromParts;
2347  delete []underloadedPartIndices;
2348  delete []overloadedPartIndices;
2349  delete []extraInPart;
2350  delete []partNext;
2351  delete []partBegins;
2352  delete []numPointsInPartsInclusiveUptoMyIndex;
2353  delete []numPointsInParts;
2354  delete []numGlobalPointsInParts;
2355 
2356  delete []shiftProcCount;
2357  delete []bestDimProcs;
2358  delete []dim_slices;
2359  delete []minCoords;
2360  delete []maxCoords;
2361  }
2362 
2363  //given the part numbers for each local coordinate,
2364  //distributes the coordinates to the corresponding processors.
2365  void distribute_points(int *coordinate_grid_parts){
2366 
2367  Tpetra::Distributor distributor(comm);
2368  ArrayView<const int> pIds( coordinate_grid_parts, this->numLocalCoords);
2369  /*
2370  for (int i = 0 ; i < this->numLocalCoords; ++i){
2371  cout << "me:" << this->myRank << " to part:" << coordinate_grid_parts[i] << endl;
2372  }
2373  */
2374  gno_t numMyNewGnos = distributor.createFromSends(pIds);
2375 
2376  //cout << "distribution step 1 me:" << this->myRank << " numLocal:" <<numMyNewGnos << " old:" << numLocalCoords << endl;
2377 
2378  this->numLocalCoords = numMyNewGnos;
2379 
2380 
2381  ArrayRCP<scalar_t> recvBuf2(distributor.getTotalReceiveLength());
2382 
2383  for (int i = 0; i < this->coordinate_dimension; ++i){
2384  ArrayView<scalar_t> s(this->coords[i], this->numLocalCoords);
2385  distributor.doPostsAndWaits<scalar_t>(s, 1, recvBuf2());
2386  delete [] this->coords[i];
2387  this->coords[i] = new scalar_t[this->numLocalCoords];
2388  for (lno_t j = 0; j < this->numLocalCoords; ++j){
2389  this->coords[i][j] = recvBuf2[j];
2390  }
2391 
2392  }
2393  }
2394 
2395  //calls MJ for p = numProcs
2396  int predistributeMJ(int *coordinate_grid_parts){
2397  int coord_dim = this->coordinate_dimension;
2398 
2399  lno_t numLocalPoints = this->numLocalCoords;
2400  gno_t numGlobalPoints = this->numGlobalCoords;
2401 
2402 
2403  //T **weight = NULL;
2404  //typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> tMVector_t;
2405  RCP<Tpetra::Map<lno_t, gno_t, node_t> > mp = rcp(
2406  new Tpetra::Map<lno_t, gno_t, node_t> (numGlobalPoints, numLocalPoints, 0, comm));
2407 
2408  Teuchos::Array<Teuchos::ArrayView<const scalar_t> > coordView(coord_dim);
2409 
2410 
2411 
2412  for (int i=0; i < coord_dim; i++){
2413  if(numLocalPoints > 0){
2414  Teuchos::ArrayView<const scalar_t> a(coords[i], numLocalPoints);
2415  coordView[i] = a;
2416  } else{
2417  Teuchos::ArrayView<const scalar_t> a;
2418  coordView[i] = a;
2419  }
2420  }
2421 
2422  RCP< Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >tmVector = RCP< Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >(
2423  new Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t>( mp, coordView.view(0, coord_dim), coord_dim));
2424 
2425 
2426  RCP<const tMVector_t> coordsConst = Teuchos::rcp_const_cast<const tMVector_t>(tmVector);
2427  std::vector<const scalar_t *> weights;
2428  std::vector <int> stride;
2429 
2430  typedef Zoltan2::XpetraMultiVectorAdapter<tMVector_t> inputAdapter_t;
2431  //inputAdapter_t ia(coordsConst);
2432  inputAdapter_t ia(coordsConst,weights, stride);
2433 
2434  Teuchos::RCP <Teuchos::ParameterList> params ;
2435  params =RCP <Teuchos::ParameterList> (new Teuchos::ParameterList, true);
2436 
2437 
2438  params->set("algorithm", "multijagged");
2439  params->set("num_global_parts", this->worldSize);
2440 
2441  //TODO we need to fix the setting parts.
2442  //Although MJ sets the parts with
2443  //currently the part setting is not correct when migration is done.
2444  params->set("migration_check_option", 2);
2445 
2446 
2448 
2449 
2450  try {
2451 #ifdef HAVE_ZOLTAN2_MPI
2452  problem = new Zoltan2::PartitioningProblem<inputAdapter_t>(&ia, params.getRawPtr(),
2453  MPI_COMM_WORLD);
2454 #else
2455  problem = new Zoltan2::PartitioningProblem<inputAdapter_t>(&ia, params.getRawPtr());
2456 #endif
2457  }
2458  CATCH_EXCEPTIONS("PartitioningProblem()")
2459 
2460  try {
2461  problem->solve();
2462  }
2463  CATCH_EXCEPTIONS("solve()")
2464 
2465  const typename inputAdapter_t::part_t *partIds = problem->getSolution().getPartListView();
2466 
2467  for (lno_t i = 0; i < this->numLocalCoords;++i){
2468  coordinate_grid_parts[i] = partIds[i];
2469  //cout << "me:" << this->myRank << " i:" << i << " goes to:" << partIds[i] << endl;
2470  }
2471  delete problem;
2472  return 0;
2473  }
2474 
2475 #ifdef HAVE_ZOLTAN2_ZOLTAN
2476  //calls RCP for p = numProcs
2477  int predistributeRCB(int *coordinate_grid_parts){
2478  int rank = this->myRank;
2479  int nprocs = this->worldSize;
2480  DOTS<tMVector_t> dots_;
2481 
2482  MEMORY_CHECK(rank==0 || rank==nprocs-1, "After initializing MPI");
2483 
2484 
2485  int nWeights = 0;
2486  int debugLevel=0;
2487  string memoryOn("memoryOn");
2488  string memoryOff("memoryOff");
2489  bool doMemory=false;
2490  int numGlobalParts = nprocs;
2491  int dummyTimer=0;
2492  bool remap=0;
2493 
2494  string balanceCount("balance_object_count");
2495  string balanceWeight("balance_object_weight");
2496  string mcnorm1("multicriteria_minimize_total_weight");
2497  string mcnorm2("multicriteria_balance_total_maximum");
2498  string mcnorm3("multicriteria_minimize_maximum_weight");
2499  string objective(balanceWeight); // default
2500 
2501  // Process command line input
2502  CommandLineProcessor commandLine(false, true);
2503  //commandLine.setOption("size", &numGlobalCoords,
2504  // "Approximate number of global coordinates.");
2505  int input_option = 0;
2506  commandLine.setOption("input_option", &input_option,
2507  "whether to use mesh creation, geometric generator, or file input");
2508  string inputFile = "";
2509 
2510  commandLine.setOption("input_file", &inputFile,
2511  "the input file for geometric generator or file input");
2512 
2513 
2514  commandLine.setOption("size", &numGlobalCoords,
2515  "Approximate number of global coordinates.");
2516  commandLine.setOption("numParts", &numGlobalParts,
2517  "Number of parts (default is one per proc).");
2518  commandLine.setOption("nWeights", &nWeights,
2519  "Number of weights per coordinate, zero implies uniform weights.");
2520  commandLine.setOption("debug", &debugLevel, "Zoltan1 debug level");
2521  commandLine.setOption("remap", "no-remap", &remap,
2522  "Zoltan1 REMAP parameter; disabled by default for scalability testing");
2523  commandLine.setOption("timers", &dummyTimer, "ignored");
2524  commandLine.setOption(memoryOn.c_str(), memoryOff.c_str(), &doMemory,
2525  "do memory profiling");
2526 
2527  string doc(balanceCount);
2528  doc.append(": ignore weights\n");
2529  doc.append(balanceWeight);
2530  doc.append(": balance on first weight\n");
2531  doc.append(mcnorm1);
2532  doc.append(": given multiple weights, balance their total.\n");
2533  doc.append(mcnorm3);
2534  doc.append(": given multiple weights, "
2535  "balance the maximum for each coordinate.\n");
2536  doc.append(mcnorm2);
2537  doc.append(": given multiple weights, balance the L2 norm of the weights.\n");
2538  commandLine.setOption("objective", &objective, doc.c_str());
2539 
2540  CommandLineProcessor::EParseCommandLineReturn rc =
2541  commandLine.parse(0, NULL);
2542 
2543 
2544 
2545  if (rc != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
2546  if (rc == Teuchos::CommandLineProcessor::PARSE_HELP_PRINTED) {
2547  if (rank==0) cout << "PASS" << endl;
2548  return 1;
2549  }
2550  else {
2551  if (rank==0) cout << "FAIL" << endl;
2552  return 0;
2553  }
2554  }
2555 
2556  //MEMORY_CHECK(doMemory && rank==0, "After processing parameters");
2557 
2558  // Create the data structure
2559 
2560  int coord_dim = this->coordinate_dimension;
2561 
2562 
2563  RCP<Tpetra::Map<lno_t, gno_t, node_t> > mp = rcp(
2564  new Tpetra::Map<lno_t, gno_t, node_t> (this->numGlobalCoords, this->numLocalCoords, 0, this->comm));
2565 
2566  Teuchos::Array<Teuchos::ArrayView<const scalar_t> > coordView(coord_dim);
2567  for (int i=0; i < coord_dim; i++){
2568  if(numLocalCoords > 0){
2569  Teuchos::ArrayView<const scalar_t> a(coords[i], numLocalCoords);
2570  coordView[i] = a;
2571  } else{
2572  Teuchos::ArrayView<const scalar_t> a;
2573  coordView[i] = a;
2574  }
2575  }
2576 
2577  tMVector_t *tmVector = new tMVector_t( mp, coordView.view(0, coord_dim), coord_dim);
2578 
2579  dots_.coordinates = tmVector;
2580  dots_.weights.resize(nWeights);
2581 
2582 
2583  MEMORY_CHECK(doMemory && rank==0, "After creating input");
2584 
2585  // Now call Zoltan to partition the problem.
2586 
2587  float ver;
2588  int aok = Zoltan_Initialize(0,NULL, &ver);
2589 
2590  if (aok != 0){
2591  printf("Zoltan_Initialize failed\n");
2592  exit(0);
2593  }
2594 
2595  struct Zoltan_Struct *zz;
2596  zz = Zoltan_Create(MPI_COMM_WORLD);
2597 
2598  Zoltan_Set_Param(zz, "LB_METHOD", "RCB");
2599  Zoltan_Set_Param(zz, "LB_APPROACH", "PARTITION");
2600  Zoltan_Set_Param(zz, "CHECK_GEOM", "0");
2601  Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1");
2602  Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "0");
2603  Zoltan_Set_Param(zz, "RETURN_LISTS", "PART");
2604  std::ostringstream oss;
2605  oss << numGlobalParts;
2606  Zoltan_Set_Param(zz, "NUM_GLOBAL_PARTS", oss.str().c_str());
2607  oss.str("");
2608  oss << debugLevel;
2609  Zoltan_Set_Param(zz, "DEBUG_LEVEL", oss.str().c_str());
2610 
2611  if (remap)
2612  Zoltan_Set_Param(zz, "REMAP", "1");
2613  else
2614  Zoltan_Set_Param(zz, "REMAP", "0");
2615 
2616  if (objective != balanceCount){
2617  oss.str("");
2618  oss << nWeights;
2619  Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", oss.str().c_str());
2620 
2621  if (objective == mcnorm1)
2622  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "1");
2623  else if (objective == mcnorm2)
2624  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "2");
2625  else if (objective == mcnorm3)
2626  Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "3");
2627  }
2628  else{
2629  Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", "0");
2630  }
2631 
2632  Zoltan_Set_Num_Obj_Fn(zz, getNumObj<tMVector_t>, &dots_);
2633  Zoltan_Set_Obj_List_Fn(zz, getObjList<tMVector_t>, &dots_);
2634  Zoltan_Set_Num_Geom_Fn(zz, getDim<tMVector_t>, &dots_);
2635  Zoltan_Set_Geom_Multi_Fn(zz, getCoords<tMVector_t>, &dots_);
2636 
2637  int changes, numGidEntries, numLidEntries, numImport, numExport;
2638  ZOLTAN_ID_PTR importGlobalGids, importLocalGids;
2639  ZOLTAN_ID_PTR exportGlobalGids, exportLocalGids;
2640  int *importProcs, *importToPart, *exportProcs, *exportToPart;
2641 
2642  MEMORY_CHECK(doMemory && rank==0, "Before Zoltan_LB_Partition");
2643 
2644 
2645  aok = Zoltan_LB_Partition(zz, &changes, &numGidEntries, &numLidEntries,
2646  &numImport, &importGlobalGids, &importLocalGids,
2647  &importProcs, &importToPart,
2648  &numExport, &exportGlobalGids, &exportLocalGids,
2649  &exportProcs, &exportToPart);
2650 
2651 
2652  MEMORY_CHECK(doMemory && rank==0, "After Zoltan_LB_Partition");
2653 
2654  for (lno_t i = 0; i < numLocalCoords; i++)
2655  coordinate_grid_parts[i] = exportToPart[i];
2656  Zoltan_Destroy(&zz);
2657  MEMORY_CHECK(doMemory && rank==0, "After Zoltan_Destroy");
2658 
2659  delete dots_.coordinates;
2660  return 0;
2661 }
2662 #endif
2664  int *coordinate_grid_parts = new int[this->numLocalCoords];
2665  switch (this->predistribution){
2666  case 1:
2667 #ifdef HAVE_ZOLTAN2_ZOLTAN
2668  this->predistributeRCB(coordinate_grid_parts);
2669  break;
2670 #endif
2671  case 2:
2672 
2673  this->predistributeMJ(coordinate_grid_parts);
2674  break;
2675  case 3:
2676  //block
2677  blockPartition(coordinate_grid_parts);
2678  break;
2679  }
2680  this->distribute_points(coordinate_grid_parts);
2681 
2682  delete []coordinate_grid_parts;
2683 
2684 
2685  }
2686 
2687  //############################################################//
2689  //############################################################//
2690 
2691 
2693  return this->numWeightsPerCoord;
2694  }
2696  return this->coordinate_dimension;
2697  }
2699  return this->numLocalCoords;
2700  }
2702  return this->numGlobalCoords;
2703  }
2704 
2706  return this->coords;
2707  }
2708 
2709  scalar_t **getLocalWeightsView(){
2710  return this->wghts;
2711  }
2712 
2713  void getLocalCoordinatesCopy( scalar_t ** c){
2714  for(int ii = 0; ii < this->coordinate_dimension; ++ii){
2715 #ifdef HAVE_ZOLTAN2_OMP
2716 #pragma omp parallel for
2717 #endif
2718  for (lno_t i = 0; i < this->numLocalCoords; ++i){
2719  c[ii][i] = this->coords[ii][i];
2720  }
2721  }
2722  }
2723 
2724  void getLocalWeightsCopy(scalar_t **w){
2725  for(int ii = 0; ii < this->numWeightsPerCoord; ++ii){
2726 #ifdef HAVE_ZOLTAN2_OMP
2727 #pragma omp parallel for
2728 #endif
2729  for (lno_t i = 0; i < this->numLocalCoords; ++i){
2730  w[ii][i] = this->wghts[ii][i];
2731  }
2732  }
2733  }
2734 };
2735 }
2736 
2737 #endif /* GEOMETRICGENERATOR */
RectangleHole(CoordinatePoint< T > center_, T edge_x_, T edge_y_)
virtual weighttype get1DWeight(T x)
virtual bool isInArea(CoordinatePoint< T > dot)
#define INVALID(STR)
CoordinateGridDistribution(gno_t alongX, gno_t alongY, gno_t alongZ, int dim, T l_x, T r_x, T l_y, T r_y, T l_z, T r_z, int myRank_, int wSize)
virtual weighttype get1DWeight(T x)=0
paramName
Definition: xml2dox.py:207
void GetPoints(lno_t requestedPointcount, CoordinatePoint< T > *points, Hole< T > **holes, lno_t holeCount, float *sharedRatios_, int myRank)
virtual bool isInArea(CoordinatePoint< T > dot)
virtual bool isInArea(CoordinatePoint< T > dot)
virtual CoordinatePoint< T > getPoint(gno_t pindex, unsigned int &state)
virtual weighttype get2DWeight(T x, T y)
void GetPoints(lno_t requestedPointcount, T **coords, lno_t tindex, Hole< T > **holes, lno_t holeCount, float *sharedRatios_, int myRank)
Hole(CoordinatePoint< T > center_, T edge1_, T edge2_, T edge3_)
CoordinateUniformDistribution(gno_t np_, int dim, T l_x, T r_x, T l_y, T r_y, T l_z, T r_z, int wSize)
SphereHole(CoordinatePoint< T > center_, T edge_)
CircleHole(CoordinatePoint< T > center_, T edge_)
Defines the PartitioningSolution class.
void perturb_data(float used_perturbation_ratio, int scale)
Start perturbation function##########################//
virtual weighttype get3DWeight(T x, T y, T z)
virtual weighttype get3DWeight(T x, T y, T z)=0
Defines the XpetraMultiVectorAdapter.
void blockPartition(int *coordinate_grid_parts)
RectangularPrismHole(CoordinatePoint< T > center_, T edge_x_, T edge_y_, T edge_z_)
CoordinateNormalDistribution(gno_t np_, int dim, CoordinatePoint< T > center_, T sd_x, T sd_y, T sd_z, int wSize)
virtual bool isInArea(CoordinatePoint< T > dot)
virtual bool isInArea(CoordinatePoint< T > dot)
SteppedEquation(T a1_, T a2_, T a3_, T b1_, T b2_, T b3_, T c_, T x1_, T y1_, T z1_, T *steps_, T *values_, int stepCount_)
const std::string shapes[]
virtual weighttype getWeight(CoordinatePoint< T > p)
outfile
Definition: xml2dox.py:11
const std::string weight_distribution_string
virtual CoordinatePoint< T > getPoint(gno_t pindex, unsigned int &state)
virtual weighttype get2DWeight(T x, T y)=0
An adapter for Xpetra::MultiVector.
void getBestSurface(int remaining, int *dimProcs, int dim, int currentDim, int &bestSurface, int *bestDimProcs)
Start Predistribution functions######################//
Expression is a generic following method.
std::string toString(tt obj)
Converts the given object to string.
void distribute_points(int *coordinate_grid_parts)
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
#define MAX_ITER_ALLOWED
GeometricGenerator(Teuchos::ParameterList &params, const RCP< const Teuchos::Comm< int > > &comm_)
Tpetra::MultiVector< double, int, int > tMVector_t
int getNumWeights()
##END Predistribution functions######################//
CoordinatePoint< T > center
#define INVALID_SHAPE_ARG(SHAPE, REQUIRED)
#define INVALIDSHAPE(STR, DIM)
#define CATCH_EXCEPTIONS(pp)
Defines the PartitioningProblem class.
CoordinateDistribution(gno_t np_, int dim, int wSize)
int predistributeMJ(int *coordinate_grid_parts)
#define epsilon
SquareHole(CoordinatePoint< T > center_, T edge_)
void solve(bool updateInputData=true)
Direct the problem to create a solution.
virtual CoordinatePoint< T > getPoint(gno_t pindex, unsigned int &state)
CubeHole(CoordinatePoint< T > center_, T edge_)
#define MAX_WEIGHT_DIM
virtual bool isInArea(CoordinatePoint< T > dot)
void getMinMaxCoords(scalar_t *globalMaxCoords, scalar_t *globalMinCoords)
#define MEMORY_CHECK(iPrint, msg)