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