Zoltan2
MultiJaggedTest.cpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
51 #include <Zoltan2_TestHelpers.hpp>
56 #include <GeometricGenerator.hpp>
57 #include <vector>
58 
60 
61 #include "Teuchos_XMLParameterListHelpers.hpp"
62 
63 #include <Teuchos_LAPACK.hpp>
64 #include <fstream>
65 #include <string>
66 using namespace std;
67 using Teuchos::RCP;
68 using Teuchos::rcp;
69 
70 
71 //#define hopper_separate_test
72 #ifdef hopper_separate_test
73 #include "stdio.h"
74 #endif
75 #define CATCH_EXCEPTIONS_AND_RETURN(pp) \
76  catch (std::runtime_error &e) { \
77  cout << "Runtime exception returned from " << pp << ": " \
78  << e.what() << " FAIL" << endl; \
79  return -1; \
80  } \
81  catch (std::logic_error &e) { \
82  cout << "Logic exception returned from " << pp << ": " \
83  << e.what() << " FAIL" << endl; \
84  return -1; \
85  } \
86  catch (std::bad_alloc &e) { \
87  cout << "Bad_alloc exception returned from " << pp << ": " \
88  << e.what() << " FAIL" << endl; \
89  return -1; \
90  } \
91  catch (std::exception &e) { \
92  cout << "Unknown exception returned from " << pp << ": " \
93  << e.what() << " FAIL" << endl; \
94  return -1; \
95  }
96 
97 #define CATCH_EXCEPTIONS_WITH_COUNT(ierr, pp) \
98  catch (std::runtime_error &e) { \
99  cout << "Runtime exception returned from " << pp << ": " \
100  << e.what() << " FAIL" << endl; \
101  (ierr)++; \
102  } \
103  catch (std::logic_error &e) { \
104  cout << "Logic exception returned from " << pp << ": " \
105  << e.what() << " FAIL" << endl; \
106  (ierr)++; \
107  } \
108  catch (std::bad_alloc &e) { \
109  cout << "Bad_alloc exception returned from " << pp << ": " \
110  << e.what() << " FAIL" << endl; \
111  (ierr)++; \
112  } \
113  catch (std::exception &e) { \
114  cout << "Unknown exception returned from " << pp << ": " \
115  << e.what() << " FAIL" << endl; \
116  (ierr)++; \
117  }
118 
119 
120 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> tMVector_t;
121 
127 const char param_comment = '#';
128 
130  const string& s,
131  const string& delimiters = " \f\n\r\t\v" )
132 {
133  return s.substr( 0, s.find_last_not_of( delimiters ) + 1 );
134 }
135 
137  const string& s,
138  const string& delimiters = " \f\n\r\t\v" )
139 {
140  return s.substr( s.find_first_not_of( delimiters ) );
141 }
142 
143 string trim_copy(
144  const string& s,
145  const string& delimiters = " \f\n\r\t\v" )
146 {
147  return trim_left_copy( trim_right_copy( s, delimiters ), delimiters );
148 }
149 
150 template <typename Adapter>
152  const char *str,
153  int dim,
154  typename Adapter::scalar_t *lower,
155  typename Adapter::scalar_t *upper,
156  size_t nparts,
157  typename Adapter::part_t *parts
158 )
159 {
160  std::cout << "boxAssign test " << str << ": Box (";
161  for (int j = 0; j < dim; j++) std::cout << lower[j] << " ";
162  std::cout << ") x (";
163  for (int j = 0; j < dim; j++) std::cout << upper[j] << " ";
164 
165  if (nparts == 0)
166  std::cout << ") does not overlap any parts" << std::endl;
167  else {
168  std::cout << ") overlaps parts ";
169  for (size_t k = 0; k < nparts; k++) std::cout << parts[k] << " ";
170  std::cout << std::endl;
171  }
172 }
173 
174 template <typename Adapter>
177  RCP<tMVector_t> &coords)
178 {
179  int ierr = 0;
180 
181  // pointAssign tests
182  int coordDim = coords->getNumVectors();
183  zscalar_t *pointDrop = new zscalar_t[coordDim];
184  typename Adapter::part_t part = -1;
185 
186  char mechar[10];
187  sprintf(mechar, "%d", problem->getComm()->getRank());
188  string me(mechar);
189 
190  // test correctness of pointAssign for owned points
191  {
192  const typename Adapter::part_t *solnPartView =
193  problem->getSolution().getPartListView();
194 
195  size_t numPoints = coords->getLocalLength();
196  for (size_t localID = 0; localID < numPoints; localID++) {
197 
198  typename Adapter::part_t solnPart = solnPartView[localID];
199 
200  for (int i = 0; i < coordDim; i++)
201  pointDrop[i] = coords->getData(i)[localID];
202 
203  try {
204  part = problem->getSolution().pointAssign(coordDim, pointDrop);
205  }
206  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + ": pointAssign -- OwnedPoints");
207 
208  std::cout << me << " Point " << localID
209  << " gid " << coords->getMap()->getGlobalElement(localID)
210  << " (" << pointDrop[0];
211  if (coordDim > 1) std::cout << " " << pointDrop[1];
212  if (coordDim > 2) std::cout << " " << pointDrop[2];
213  std::cout << ") in boxPart " << part
214  << " in solnPart " << solnPart
215  << std::endl;
216 
217 // this error test does not work for points that fall on the cuts.
218 // like Zoltan's RCB, pointAssign arbitrarily picks a part along the cut.
219 // the arbitrarily chosen part will not necessarily be the one to which
220 // the coordinate was assigned in partitioning.
221 //
222 // if (part != solnPart) {
223 // std::cout << me << " pointAssign: incorrect part " << part
224 // << " found; should be " << solnPart
225 // << " for point " << j << std::endl;
226 // ierr++;
227 // }
228  }
229  }
230 
231  {
233  typename Adapter::part_t> >
234  pBoxes = problem->getSolution().getPartBoxesView();
235  for (size_t i = 0; i < pBoxes.size(); i++) {
236  zscalar_t *lmin = pBoxes[i].getlmins();
237  zscalar_t *lmax = pBoxes[i].getlmaxs();;
238  std::cout << me << " pBox " << i << " pid " << pBoxes[i].getpId()
239  << " (" << lmin[0] << "," << lmin[1] << ","
240  << (coordDim > 2 ? lmin[2] : 0) << ") x "
241  << " (" << lmax[0] << "," << lmax[1] << ","
242  << (coordDim > 2 ? lmax[2] : 0) << ")" << std::endl;
243  }
244  }
245 
246  // test the origin
247  {
248  for (int i = 0; i < coordDim; i++) pointDrop[i] = 0.;
249  try {
250  part = problem->getSolution().pointAssign(coordDim, pointDrop);
251  }
252  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " pointAssign -- Origin");
253  std::cout << me << " OriginPoint (" << pointDrop[0];
254  if (coordDim > 1) std::cout << " " << pointDrop[1];
255  if (coordDim > 2) std::cout << " " << pointDrop[2];
256  std::cout << ") part " << part << std::endl;
257  }
258 
259  // test point with negative coordinates
260  {
261  for (int i = 0; i < coordDim; i++) pointDrop[i] = -100.+i;
262  try {
263  part = problem->getSolution().pointAssign(coordDim, pointDrop);
264  }
265  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " pointAssign -- Negative Point");
266  std::cout << me << " NegativePoint (" << pointDrop[0];
267  if (coordDim > 1) std::cout << " " << pointDrop[1];
268  if (coordDim > 2) std::cout << " " << pointDrop[2];
269  std::cout << ") part " << part << std::endl;
270  }
271 
272  // test a point that's way out there
273  {
274  for (int i = 0; i < coordDim; i++) pointDrop[i] = i*5;
275  try {
276  part = problem->getSolution().pointAssign(coordDim, pointDrop);
277  }
278  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " pointAssign -- i*5 Point");
279  std::cout << me << " i*5-Point (" << pointDrop[0];
280  if (coordDim > 1) std::cout << " " << pointDrop[1];
281  if (coordDim > 2) std::cout << " " << pointDrop[2];
282  std::cout << ") part " << part << std::endl;
283  }
284 
285  // test a point that's way out there
286  {
287  for (int i = 0; i < coordDim; i++) pointDrop[i] = 10+i*5;
288  try {
289  part = problem->getSolution().pointAssign(coordDim, pointDrop);
290  }
291  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " pointAssign -- WoopWoop");
292  std::cout << me << " WoopWoop-Point (" << pointDrop[0];
293  if (coordDim > 1) std::cout << " " << pointDrop[1];
294  if (coordDim > 2) std::cout << " " << pointDrop[2];
295  std::cout << ") part " << part << std::endl;
296  }
297 
298  delete [] pointDrop;
299  return ierr;
300 }
301 
302 template <typename Adapter>
305  RCP<tMVector_t> &coords)
306 {
307  int ierr = 0;
308 
309  // boxAssign tests
310  int coordDim = coords->getNumVectors();
311  zscalar_t *lower = new zscalar_t[coordDim];
312  zscalar_t *upper = new zscalar_t[coordDim];
313 
314  char mechar[10];
315  sprintf(mechar, "%d", problem->getComm()->getRank());
316  string me(mechar);
317 
319  typename Adapter::part_t> >
320  pBoxes = problem->getSolution().getPartBoxesView();
321  size_t nBoxes = pBoxes.size();
322 
323  // test a box that is smaller than a part
324  {
325  size_t nparts;
326  typename Adapter::part_t *parts;
327  size_t pickabox = nBoxes / 2;
328  for (int i = 0; i < coordDim; i++) {
329  zscalar_t dd = 0.2 * (pBoxes[pickabox].getlmaxs()[i] -
330  pBoxes[pickabox].getlmins()[i]);
331  lower[i] = pBoxes[pickabox].getlmins()[i] + dd;
332  upper[i] = pBoxes[pickabox].getlmaxs()[i] - dd;
333  }
334  try {
335  problem->getSolution().boxAssign(coordDim, lower, upper,
336  nparts, &parts);
337  }
338  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " boxAssign -- smaller");
339  if (nparts > 1) {
340  std::cout << me << " FAIL boxAssign error: smaller test, nparts > 1"
341  << std::endl;
342  ierr++;
343  }
344  print_boxAssign_result<Adapter>("smallerbox", coordDim,
345  lower, upper, nparts, parts);
346  delete [] parts;
347  }
348 
349  // test a box that is larger than a part
350  {
351  size_t nparts;
352  typename Adapter::part_t *parts;
353  size_t pickabox = nBoxes / 2;
354  for (int i = 0; i < coordDim; i++) {
355  zscalar_t dd = 0.2 * (pBoxes[pickabox].getlmaxs()[i] -
356  pBoxes[pickabox].getlmins()[i]);
357  lower[i] = pBoxes[pickabox].getlmins()[i] - dd;
358  upper[i] = pBoxes[pickabox].getlmaxs()[i] + dd;
359  }
360  try {
361  problem->getSolution().boxAssign(coordDim, lower, upper,
362  nparts, &parts);
363  }
364  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " boxAssign -- larger");
365 
366  // larger box should have at least two parts in it for k > 1.
367  if ((nBoxes > 1) && (nparts < 2)) {
368  std::cout << me << " FAIL boxAssign error: "
369  << "larger test, nparts < 2"
370  << std::endl;
371  ierr++;
372  }
373 
374  // parts should include pickabox's part
375  bool found_pickabox = 0;
376  for (size_t i = 0; i < nparts; i++)
377  if (parts[i] == pBoxes[pickabox].getpId()) {
378  found_pickabox = 1;
379  break;
380  }
381  if (!found_pickabox) {
382  std::cout << me << " FAIL boxAssign error: "
383  << "larger test, pickabox not found"
384  << std::endl;
385  ierr++;
386  }
387 
388  print_boxAssign_result<Adapter>("largerbox", coordDim,
389  lower, upper, nparts, parts);
390  delete [] parts;
391  }
392 
393  // test a box that includes all parts
394  {
395  size_t nparts;
396  typename Adapter::part_t *parts;
397  for (int i = 0; i < coordDim; i++) {
398  lower[i] = std::numeric_limits<zscalar_t>::max();
399  upper[i] = std::numeric_limits<zscalar_t>::min();
400  }
401  for (size_t j = 0; j < nBoxes; j++) {
402  for (int i = 0; i < coordDim; i++) {
403  if (pBoxes[j].getlmins()[i] <= lower[i])
404  lower[i] = pBoxes[j].getlmins()[i];
405  if (pBoxes[j].getlmaxs()[i] >= upper[i])
406  upper[i] = pBoxes[j].getlmaxs()[i];
407  }
408  }
409  try {
410  problem->getSolution().boxAssign(coordDim, lower, upper,
411  nparts, &parts);
412  }
413  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " boxAssign -- global");
414 
415  // global box should have all parts
416  if (nparts != nBoxes) {
417  std::cout << me << " FAIL boxAssign error: "
418  << "global test, nparts found " << nparts
419  << " != num global parts " << nBoxes
420  << std::endl;
421  ierr++;
422  }
423  print_boxAssign_result<Adapter>("globalbox", coordDim,
424  lower, upper, nparts, parts);
425  delete [] parts;
426  }
427 
428  // test a box that is bigger than the entire domain
429  // Assuming lower and upper are still set to the global box boundary
430  // from the previous test
431  {
432  size_t nparts;
433  typename Adapter::part_t *parts;
434  for (int i = 0; i < coordDim; i++) {
435  lower[i] -= 2.;
436  upper[i] += 2.;
437  }
438 
439  try {
440  problem->getSolution().boxAssign(coordDim, lower, upper,
441  nparts, &parts);
442  }
443  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " boxAssign -- bigdomain");
444 
445  // bigdomain box should have all parts
446  if (nparts != nBoxes) {
447  std::cout << me << " FAIL boxAssign error: "
448  << "bigdomain test, nparts found " << nparts
449  << " != num global parts " << nBoxes
450  << std::endl;
451  ierr++;
452  }
453  print_boxAssign_result<Adapter>("bigdomainbox", coordDim,
454  lower, upper, nparts, parts);
455  delete [] parts;
456  }
457 
458  // test a box that is way out there
459  // Assuming lower and upper are still set to at least the global box
460  // boundary from the previous test
461  {
462  size_t nparts;
463  typename Adapter::part_t *parts;
464  for (int i = 0; i < coordDim; i++) {
465  lower[i] = upper[i] + 10;
466  upper[i] += 20;
467  }
468 
469  try {
470  problem->getSolution().boxAssign(coordDim, lower, upper,
471  nparts, &parts);
472  }
473  CATCH_EXCEPTIONS_WITH_COUNT(ierr, me + " boxAssign -- out there");
474 
475  // For now, boxAssign returns zero if there is no box overlap.
476  // TODO: this result should be changed in boxAssign definition
477  if (nparts != 0) {
478  std::cout << me << " FAIL boxAssign error: "
479  << "outthere test, nparts found " << nparts
480  << " != zero"
481  << std::endl;
482  ierr++;
483  }
484  print_boxAssign_result<Adapter>("outthere box", coordDim,
485  lower, upper, nparts, parts);
486  delete [] parts;
487  }
488 
489  delete [] lower;
490  delete [] upper;
491  return ierr;
492 }
493 
494 void readGeoGenParams(string paramFileName, Teuchos::ParameterList &geoparams, const RCP<const Teuchos::Comm<int> > & comm){
495  std::string input = "";
496  char inp[25000];
497  for(int i = 0; i < 25000; ++i){
498  inp[i] = 0;
499  }
500 
501  bool fail = false;
502  if(comm->getRank() == 0){
503 
504  fstream inParam(paramFileName.c_str());
505  if (inParam.fail())
506  {
507  fail = true;
508  }
509  if(!fail)
510  {
511  std::string tmp = "";
512  getline (inParam,tmp);
513  while (!inParam.eof()){
514  if(tmp != ""){
515  tmp = trim_copy(tmp);
516  if(tmp != ""){
517  input += tmp + "\n";
518  }
519  }
520  getline (inParam,tmp);
521  }
522  inParam.close();
523  for (size_t i = 0; i < input.size(); ++i){
524  inp[i] = input[i];
525  }
526  }
527  }
528 
529 
530 
531  int size = input.size();
532  if(fail){
533  size = -1;
534  }
535  comm->broadcast(0, sizeof(int), (char*) &size);
536  if(size == -1){
537  throw "File " + paramFileName + " cannot be opened.";
538  }
539  comm->broadcast(0, size, inp);
540  istringstream inParam(inp);
541  string str;
542  getline (inParam,str);
543  while (!inParam.eof()){
544  if(str[0] != param_comment){
545  size_t pos = str.find('=');
546  if(pos == string::npos){
547  throw "Invalid Line:" + str + " in parameter file";
548  }
549  string paramname = trim_copy(str.substr(0,pos));
550  string paramvalue = trim_copy(str.substr(pos + 1));
551  geoparams.set(paramname, paramvalue);
552  }
553  getline (inParam,str);
554  }
555 }
556 
557 int GeometricGenInterface(RCP<const Teuchos::Comm<int> > &comm,
558  int numParts, float imbalance,
559  std::string paramFile, std::string pqParts,
560  std::string pfname,
561  int k,
562  int migration_check_option,
563  int migration_all_to_all_type,
564  zscalar_t migration_imbalance_cut_off,
565  int migration_processor_assignment_type,
566  int migration_doMigration_type,
567  bool test_boxes,
568  bool rectilinear
569 )
570 {
571  int ierr = 0;
572  Teuchos::ParameterList geoparams("geo params");
573  readGeoGenParams(paramFile, geoparams, comm);
576  comm);
577 
578  int coord_dim = gg->getCoordinateDimension();
579  int numWeightsPerCoord = gg->getNumWeights();
580  zlno_t numLocalPoints = gg->getNumLocalCoords();
581  zgno_t numGlobalPoints = gg->getNumGlobalCoords();
582  zscalar_t **coords = new zscalar_t * [coord_dim];
583  for(int i = 0; i < coord_dim; ++i){
584  coords[i] = new zscalar_t[numLocalPoints];
585  }
586  gg->getLocalCoordinatesCopy(coords);
587  zscalar_t **weight = NULL;
588  if (numWeightsPerCoord) {
589  weight= new zscalar_t * [numWeightsPerCoord];
590  for(int i = 0; i < numWeightsPerCoord; ++i){
591  weight[i] = new zscalar_t[numLocalPoints];
592  }
593  gg->getLocalWeightsCopy(weight);
594  }
595 
596  delete gg;
597 
598  RCP<Tpetra::Map<zlno_t, zgno_t, znode_t> > mp = rcp(
599  new Tpetra::Map<zlno_t, zgno_t, znode_t>(numGlobalPoints,
600  numLocalPoints, 0, comm));
601 
602  Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(coord_dim);
603  for (int i=0; i < coord_dim; i++){
604  if(numLocalPoints > 0){
605  Teuchos::ArrayView<const zscalar_t> a(coords[i], numLocalPoints);
606  coordView[i] = a;
607  }
608  else {
609  Teuchos::ArrayView<const zscalar_t> a;
610  coordView[i] = a;
611  }
612  }
613 
614  RCP<tMVector_t> tmVector = RCP<tMVector_t>(new
615  tMVector_t(mp, coordView.view(0, coord_dim),
616  coord_dim));
617 
618  RCP<const tMVector_t> coordsConst =
619  Teuchos::rcp_const_cast<const tMVector_t>(tmVector);
620  vector<const zscalar_t *> weights;
621  if(numWeightsPerCoord){
622  for (int i = 0; i < numWeightsPerCoord;++i){
623  weights.push_back(weight[i]);
624  }
625  }
626  vector <int> stride;
627 
628  typedef Zoltan2::XpetraMultiVectorAdapter<tMVector_t> inputAdapter_t;
630  //inputAdapter_t ia(coordsConst);
631  inputAdapter_t *ia = new inputAdapter_t(coordsConst,weights, stride);
632 
633  Teuchos::RCP<Teuchos::ParameterList> params ;
634 
635  //Teuchos::ParameterList params("test params");
636  if(pfname != ""){
637  params = Teuchos::getParametersFromXmlFile(pfname);
638  }
639  else {
640  params =RCP<Teuchos::ParameterList>(new Teuchos::ParameterList, true);
641  }
642 /*
643  params->set("memory_output_stream" , "std::cout");
644  params->set("memory_procs" , 0);
645  */
646  params->set("timer_output_stream" , "std::cout");
647 
648  params->set("algorithm", "multijagged");
649  if (test_boxes)
650  params->set("mj_keep_part_boxes", true); // bool parameter
651  if (rectilinear)
652  params->set("rectilinear", true ); // bool parameter
653 
654  if(imbalance > 1)
655  params->set("imbalance_tolerance", double(imbalance));
656 
657  if(pqParts != "")
658  params->set("mj_parts", pqParts);
659  if(numParts > 0)
660  params->set("num_global_parts", numParts);
661  if (k > 0)
662  params->set("mj_concurrent_part_count", k);
663  if(migration_check_option >= 0)
664  params->set("mj_migration_option", migration_check_option);
665  if(migration_imbalance_cut_off >= 0)
666  params->set("mj_minimum_migration_imbalance",
667  double(migration_imbalance_cut_off));
668 
670  try {
672  params.getRawPtr(),
673  comm);
674  }
675  CATCH_EXCEPTIONS_AND_RETURN("PartitioningProblem()")
676 
677  try {
678  problem->solve();
679  }
680  CATCH_EXCEPTIONS_AND_RETURN("solve()")
681 
682  // create metric object
683 
684  RCP<quality_t> metricObject =
685  rcp(new quality_t(ia,params.getRawPtr(),comm,&problem->getSolution()));
686 
687  if (comm->getRank() == 0){
688  metricObject->printMetrics(cout);
689  }
690  problem->printTimers();
691 
692  // run pointAssign tests
693  if (test_boxes) {
694  ierr = run_pointAssign_tests<inputAdapter_t>(problem, tmVector);
695  ierr += run_boxAssign_tests<inputAdapter_t>(problem, tmVector);
696  }
697 
698  if(numWeightsPerCoord){
699  for(int i = 0; i < numWeightsPerCoord; ++i)
700  delete [] weight[i];
701  delete [] weight;
702  }
703  if(coord_dim){
704  for(int i = 0; i < coord_dim; ++i)
705  delete [] coords[i];
706  delete [] coords;
707  }
708  delete problem;
709  delete ia;
710  return ierr;
711 }
712 
714  RCP<const Teuchos::Comm<int> > &comm,
715  int numParts,
716  float imbalance,
717  std::string fname,
718  std::string pqParts,
719  std::string pfname,
720  int k,
721  int migration_check_option,
722  int migration_all_to_all_type,
723  zscalar_t migration_imbalance_cut_off,
724  int migration_processor_assignment_type,
725  int migration_doMigration_type,
726  bool test_boxes,
727  bool rectilinear
728 )
729 {
730  int ierr = 0;
731  //std::string fname("simple");
732  //cout << "running " << fname << endl;
733 
734  UserInputForTests uinput(testDataFilePath, fname, comm, true);
735 
736  RCP<tMVector_t> coords = uinput.getUICoordinates();
737 
738  RCP<const tMVector_t> coordsConst = rcp_const_cast<const tMVector_t>(coords);
739  typedef Zoltan2::XpetraMultiVectorAdapter<tMVector_t> inputAdapter_t;
741  inputAdapter_t *ia = new inputAdapter_t(coordsConst);
742 
743  Teuchos::RCP <Teuchos::ParameterList> params ;
744 
745  //Teuchos::ParameterList params("test params");
746  if(pfname != ""){
747  params = Teuchos::getParametersFromXmlFile(pfname);
748  }
749  else {
750  params =RCP <Teuchos::ParameterList> (new Teuchos::ParameterList, true);
751  }
752 
753  //params->set("timer_output_stream" , "std::cout");
754  if (test_boxes)
755  params->set("mj_keep_part_boxes", true); // bool parameter
756  if (rectilinear)
757  params->set("rectilinear", true); // bool parameter
758  params->set("algorithm", "multijagged");
759  if(imbalance > 1){
760  params->set("imbalance_tolerance", double(imbalance));
761  }
762 
763  if(pqParts != ""){
764  params->set("mj_parts", pqParts);
765  }
766  if(numParts > 0){
767  params->set("num_global_parts", numParts);
768  }
769  if (k > 0){
770  params->set("mj_concurrent_part_count", k);
771  }
772  if(migration_check_option >= 0){
773  params->set("mj_migration_option", migration_check_option);
774  }
775  if(migration_imbalance_cut_off >= 0){
776  params->set("mj_minimum_migration_imbalance",
777  double (migration_imbalance_cut_off));
778  }
779 
781  try {
783  params.getRawPtr(),
784  comm);
785  }
786  CATCH_EXCEPTIONS_AND_RETURN("PartitioningProblem()")
787 
788  try {
789  problem->solve();
790  }
791  CATCH_EXCEPTIONS_AND_RETURN("solve()")
792 
793  {
794  // Run a test with BasicVectorAdapter and xyzxyz format coordinates
795  const int bvme = comm->getRank();
796  const inputAdapter_t::lno_t bvlen =
797  inputAdapter_t::lno_t(coords->getLocalLength());
798  const size_t bvnvecs = coords->getNumVectors();
799  const size_t bvsize = coords->getNumVectors() * coords->getLocalLength();
800 
801  ArrayRCP<inputAdapter_t::scalar_t> *bvtpetravectors =
802  new ArrayRCP<inputAdapter_t::scalar_t>[bvnvecs];
803  for (size_t i = 0; i < bvnvecs; i++)
804  bvtpetravectors[i] = coords->getDataNonConst(i);
805 
806  int idx = 0;
807  inputAdapter_t::gno_t *bvgids = new
808  inputAdapter_t::gno_t[coords->getLocalLength()];
809  inputAdapter_t::scalar_t *bvcoordarr = new inputAdapter_t::scalar_t[bvsize];
810  for (inputAdapter_t::lno_t j = 0; j < bvlen; j++) {
811  bvgids[j] = coords->getMap()->getGlobalElement(j);
812  for (size_t i = 0; i < bvnvecs; i++) {
813  bvcoordarr[idx++] = bvtpetravectors[i][j];
814  }
815  }
816 
817  typedef Zoltan2::BasicUserTypes<inputAdapter_t::scalar_t,
818  inputAdapter_t::lno_t,
819  inputAdapter_t::gno_t> bvtypes_t;
820  typedef Zoltan2::BasicVectorAdapter<bvtypes_t> bvadapter_t;
821  std::vector<const inputAdapter_t::scalar_t *> bvcoords(bvnvecs);
822  std::vector<int> bvstrides(bvnvecs);
823  for (size_t i = 0; i < bvnvecs; i++) {
824  bvcoords[i] = &bvcoordarr[i];
825  bvstrides[i] = bvnvecs;
826  }
827  std::vector<const inputAdapter_t::scalar_t *> bvwgts;
828  std::vector<int> bvwgtstrides;
829 
830  bvadapter_t bvia(bvlen, bvgids, bvcoords, bvstrides,
831  bvwgts, bvwgtstrides);
832 
834  try {
835  bvproblem = new Zoltan2::PartitioningProblem<bvadapter_t>(&bvia,
836  params.getRawPtr(),
837  comm);
838  }
839  CATCH_EXCEPTIONS_AND_RETURN("PartitioningProblem()")
840 
841  try {
842  bvproblem->solve();
843  }
844  CATCH_EXCEPTIONS_AND_RETURN("solve()")
845 
846  // Compare with MultiVectorAdapter result
847  for (inputAdapter_t::lno_t i = 0; i < bvlen; i++) {
848  if (problem->getSolution().getPartListView()[i] !=
849  bvproblem->getSolution().getPartListView()[i])
850  cout << bvme << " " << i << " "
851  << coords->getMap()->getGlobalElement(i) << " " << bvgids[i]
852  << ": XMV " << problem->getSolution().getPartListView()[i]
853  << "; BMV " << bvproblem->getSolution().getPartListView()[i]
854  << " : FAIL" << endl;
855  }
856 
857  delete [] bvgids;
858  delete [] bvcoordarr;
859  delete [] bvtpetravectors;
860  delete bvproblem;
861  }
862 
863  if (coordsConst->getGlobalLength() < 40) {
864  int len = coordsConst->getLocalLength();
865  const inputAdapter_t::part_t *zparts =
866  problem->getSolution().getPartListView();
867  for (int i = 0; i < len; i++)
868  cout << comm->getRank()
869  << " lid " << i
870  << " gid " << coords->getMap()->getGlobalElement(i)
871  << " part " << zparts[i] << endl;
872  }
873 
874  // create metric object
875 
876  RCP<quality_t> metricObject =
877  rcp(new quality_t(ia,params.getRawPtr(),comm,&problem->getSolution()));
878 
879  if (comm->getRank() == 0){
880  metricObject->printMetrics(cout);
881  cout << "testFromDataFile is done " << endl;
882  }
883 
884  problem->printTimers();
885 
886  // run pointAssign tests
887  if (test_boxes) {
888  ierr = run_pointAssign_tests<inputAdapter_t>(problem, coords);
889  ierr += run_boxAssign_tests<inputAdapter_t>(problem, coords);
890  }
891 
892  delete problem;
893  delete ia;
894  return ierr;
895 }
896 
897 #ifdef hopper_separate_test
898 
899 template <typename zscalar_t, typename zlno_t>
900 void getCoords(zscalar_t **&coords, zlno_t &numLocal, int &dim, string fileName){
901  FILE *f = fopen(fileName.c_str(), "r");
902  if (f == NULL){
903  cout << fileName << " cannot be opened" << endl;
904  exit(1);
905  }
906  fscanf(f, "%d", &numLocal);
907  fscanf(f, "%d", &dim);
908  coords = new zscalar_t *[ dim];
909  for (int i = 0; i < dim; ++i){
910  coords[i] = new zscalar_t[numLocal];
911  }
912  for (int i = 0; i < dim; ++i){
913  for (zlno_t j = 0; j < numLocal; ++j){
914  fscanf(f, "%lf", &(coords[i][j]));
915  }
916  }
917  fclose(f);
918 }
919 
920 int testFromSeparateDataFiles(
921  RCP<const Teuchos::Comm<int> > &comm,
922  int numParts,
923  float imbalance,
924  std::string fname,
925  std::string pqParts,
926  std::string pfname,
927  int k,
928  int migration_check_option,
929  int migration_all_to_all_type,
930  zscalar_t migration_imbalance_cut_off,
931  int migration_processor_assignment_type,
932  int migration_doMigration_type,
933  int test_boxes,
934  bool rectilinear
935 )
936 {
937  //std::string fname("simple");
938  //cout << "running " << fname << endl;
939 
940  int ierr = 0;
941  int mR = comm->getRank();
942  if (mR == 0) cout << "size of zscalar_t:" << sizeof(zscalar_t) << endl;
943  string tFile = fname +"_" + Teuchos::toString<int>(mR) + ".mtx";
944  zscalar_t **double_coords;
945  zlno_t numLocal = 0;
946  int dim = 0;
947  getCoords<zscalar_t, zlno_t>(double_coords, numLocal, dim, tFile);
948  //UserInputForTests uinput(testDataFilePath, fname, comm, true);
949  Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(dim);
950  for (int i=0; i < dim; i++){
951  if(numLocal > 0){
952  Teuchos::ArrayView<const zscalar_t> a(double_coords[i], numLocal);
953  coordView[i] = a;
954  } else{
955  Teuchos::ArrayView<const zscalar_t> a;
956  coordView[i] = a;
957  }
958  }
959 
960  zgno_t numGlobal;
961  zgno_t nL = numLocal;
962  Teuchos::Comm<int> *tcomm = (Teuchos::Comm<int> *)comm.getRawPtr();
963 
964  reduceAll<int, zgno_t>(
965  *tcomm,
966  Teuchos::REDUCE_SUM,
967  1,
968  &nL,
969  &numGlobal
970  );
971 
972 
973  RCP<Tpetra::Map<zlno_t, zgno_t, znode_t> > mp = rcp(
974  new Tpetra::Map<zlno_t, zgno_t, znode_t> (numGlobal, numLocal, 0, comm));
975  RCP< Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> >coords = RCP< Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> >(
976  new Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>( mp, coordView.view(0, dim), dim));
977 
978 
979 
980  RCP<const tMVector_t> coordsConst = rcp_const_cast<const tMVector_t>(coords);
981 
982  typedef Zoltan2::XpetraMultiVectorInput<tMVector_t> inputAdapter_t;
984  inputAdapter_t *ia = new inputAdapter_t(coordsConst);
985 
986  Teuchos::RCP <Teuchos::ParameterList> params ;
987 
988  //Teuchos::ParameterList params("test params");
989  if(pfname != ""){
990  params = Teuchos::getParametersFromXmlFile(pfname);
991  }
992  else {
993  params =RCP <Teuchos::ParameterList> (new Teuchos::ParameterList, true);
994  }
995 
996  //params->set("timer_output_stream" , "std::cout");
997  params->set("algorithm", "multijagged");
998  if(imbalance > 1){
999  params->set("imbalance_tolerance", double(imbalance));
1000  }
1001 
1002  if(pqParts != ""){
1003  params->set("mj_parts", pqParts);
1004  }
1005  if(numParts > 0){
1006  params->set("num_global_parts", numParts);
1007  }
1008  if (k > 0){
1009  params->set("parallel_part_calculation_count", k);
1010  }
1011  if(migration_processor_assignment_type >= 0){
1012  params->set("migration_processor_assignment_type", migration_processor_assignment_type);
1013  }
1014  if(migration_check_option >= 0){
1015  params->set("migration_check_option", migration_check_option);
1016  }
1017  if(migration_all_to_all_type >= 0){
1018  params->set("migration_all_to_all_type", migration_all_to_all_type);
1019  }
1020  if(migration_imbalance_cut_off >= 0){
1021  params->set("migration_imbalance_cut_off",
1022  double (migration_imbalance_cut_off));
1023  }
1024  if (migration_doMigration_type >= 0){
1025  params->set("migration_doMigration_type", int (migration_doMigration_type));
1026  }
1027  if (test_boxes)
1028  params->set("mj_keep_part_boxes", true); // bool parameter
1029  if (rectilinear)
1030  params->set("rectilinear", true); // bool parameter
1031 
1033  try {
1034  problem =
1036  params.getRawPtr(),
1037  comm);
1038  }
1039  CATCH_EXCEPTIONS_AND_RETURN("PartitioningProblem()")
1040 
1041  try {
1042  problem->solve();
1043  }
1044  CATCH_EXCEPTIONS_AND_RETURN("solve()")
1045 
1046  if (coordsConst->getGlobalLength() < 40) {
1047  int len = coordsConst->getLocalLength();
1048  const inputAdapter_t::part_t *zparts =
1049  problem->getSolution().getPartListView();
1050  for (int i = 0; i < len; i++)
1051  cout << comm->getRank()
1052  << " gid " << coords->getMap()->getGlobalElement(i)
1053  << " part " << zparts[i] << endl;
1054  }
1055 
1056  //create metric object
1057 
1058  RCP<quality_t> metricObject =
1059  rcp(new quality_t(ia,params.getRawPtr(),comm,&problem->getSolution()));
1060 
1061  if (comm->getRank() == 0){
1062  metricObject->printMetrics(cout);
1063  cout << "testFromDataFile is done " << endl;
1064  }
1065 
1066  problem->printTimers();
1067 
1068  // run pointAssign tests
1069  if (test_boxes) {
1070  ierr = run_pointAssign_tests<inputAdapter_t>(problem, coords);
1071  ierr += run_boxAssign_tests<inputAdapter_t>(problem, coords);
1072  }
1073 
1074  delete problem;
1075  delete ia;
1076  return ierr;
1077 }
1078 #endif
1079 
1080 
1081 
1082 string convert_to_string(char *args){
1083  string tmp = "";
1084  for(int i = 0; args[i] != 0; i++)
1085  tmp += args[i];
1086  return tmp;
1087 }
1088 bool getArgumentValue(string &argumentid, double &argumentValue, string argumentline){
1089  stringstream stream(stringstream::in | stringstream::out);
1090  stream << argumentline;
1091  getline(stream, argumentid, '=');
1092  if (stream.eof()){
1093  return false;
1094  }
1095  stream >> argumentValue;
1096  return true;
1097 }
1098 
1100  int argc,
1101  char **argv,
1102  int &numParts,
1103  float &imbalance ,
1104  string &pqParts,
1105  int &opt,
1106  std::string &fname,
1107  std::string &pfname,
1108  int &k,
1109  int &migration_check_option,
1110  int &migration_all_to_all_type,
1111  zscalar_t &migration_imbalance_cut_off,
1112  int &migration_processor_assignment_type,
1113  int &migration_doMigration_type,
1114  bool &test_boxes,
1115  bool &rectilinear
1116 )
1117 {
1118  bool isCset = false;
1119  bool isPset = false;
1120  bool isFset = false;
1121  bool isPFset = false;
1122 
1123  for(int i = 0; i < argc; ++i){
1124  string tmp = convert_to_string(argv[i]);
1125  string identifier = "";
1126  long long int value = -1; double fval = -1;
1127  if(!getArgumentValue(identifier, fval, tmp)) continue;
1128  value = (long long int) (fval);
1129 
1130  if(identifier == "C"){
1131  if(value > 0){
1132  numParts=value;
1133  isCset = true;
1134  } else {
1135  throw "Invalid argument at " + tmp;
1136  }
1137  } else if(identifier == "P"){
1138  stringstream stream(stringstream::in | stringstream::out);
1139  stream << tmp;
1140  string ttmp;
1141  getline(stream, ttmp, '=');
1142  stream >> pqParts;
1143  isPset = true;
1144  }else if(identifier == "I"){
1145  if(fval > 0){
1146  imbalance=fval;
1147  } else {
1148  throw "Invalid argument at " + tmp;
1149  }
1150  } else if(identifier == "MI"){
1151  if(fval > 0){
1152  migration_imbalance_cut_off=fval;
1153  } else {
1154  throw "Invalid argument at " + tmp;
1155  }
1156  } else if(identifier == "MO"){
1157  if(value >=0 ){
1158  migration_check_option = value;
1159  } else {
1160  throw "Invalid argument at " + tmp;
1161  }
1162  } else if(identifier == "AT"){
1163  if(value >=0 ){
1164  migration_processor_assignment_type = value;
1165  } else {
1166  throw "Invalid argument at " + tmp;
1167  }
1168  }
1169 
1170  else if(identifier == "MT"){
1171  if(value >=0 ){
1172  migration_all_to_all_type = value;
1173  } else {
1174  throw "Invalid argument at " + tmp;
1175  }
1176  }
1177  else if(identifier == "DM"){
1178  if(value >=0 ){
1179  migration_doMigration_type = value;
1180  } else {
1181  throw "Invalid argument at " + tmp;
1182  }
1183  }
1184  else if(identifier == "F"){
1185  stringstream stream(stringstream::in | stringstream::out);
1186  stream << tmp;
1187  getline(stream, fname, '=');
1188 
1189  stream >> fname;
1190  isFset = true;
1191  }
1192  else if(identifier == "PF"){
1193  stringstream stream(stringstream::in | stringstream::out);
1194  stream << tmp;
1195  getline(stream, pfname, '=');
1196 
1197  stream >> pfname;
1198  isPFset = true;
1199  }
1200 
1201  else if(identifier == "O"){
1202  if(value >= 0 && value <= 3){
1203  opt = value;
1204  } else {
1205  throw "Invalid argument at " + tmp;
1206  }
1207  }
1208  else if(identifier == "K"){
1209  if(value >=0 ){
1210  k = value;
1211  } else {
1212  throw "Invalid argument at " + tmp;
1213  }
1214  }
1215  else if(identifier == "TB"){
1216  if(value >=0 ){
1217  test_boxes = (value == 0 ? false : true);
1218  } else {
1219  throw "Invalid argument at " + tmp;
1220  }
1221  }
1222  else if(identifier == "R"){
1223  if(value >=0 ){
1224  rectilinear = (value == 0 ? false : true);
1225  } else {
1226  throw "Invalid argument at " + tmp;
1227  }
1228  }
1229  else {
1230  throw "Invalid argument at " + tmp;
1231  }
1232 
1233  }
1234  if(!( (isCset || isPset || isPFset) && isFset)){
1235  throw "(C || P || PF) && F are mandatory arguments.";
1236  }
1237 
1238 }
1239 
1240 void print_usage(char *executable){
1241  cout << "\nUsage:" << endl;
1242  cout << executable << " arglist" << endl;
1243  cout << "arglist:" << endl;
1244  cout << "\tC=numParts: numParts > 0" << endl;
1245  cout << "\tP=MultiJaggedPart: Example: P=512,512" << endl;
1246  cout << "\tI=imbalance: Example I=1.03 (ignored for now.)" << endl;
1247  cout << "\tF=filePath: When O=0 the path of the coordinate input file, for O>1 the path to the geometric generator parameter file." << endl;
1248  cout << "\tO=input option: O=0 for reading coordinate from file, O>0 for generating coordinate from coordinate generator file. Default will run geometric generator." << endl;
1249  cout << "\tK=concurrent part calculation input: K>0." << endl;
1250  cout << "\tMI=migration_imbalance_cut_off: MI=1.35. " << endl;
1251  cout << "\tMT=migration_all_to_all_type: 0 for alltoallv, 1 for Zoltan_Comm, 2 for Zoltan2 Distributor object(Default 1)." << endl;
1252  cout << "\tMO=migration_check_option: 0 for decision on imbalance, 1 for forcing migration, >1 for avoiding migration. (Default-0)" << endl;
1253  cout << "\tAT=migration_processor_assignment_type. 0-for assigning procs with respect to proc ownment, otherwise, assignment with respect to proc closeness." << endl;
1254  cout << "Example:\n" << executable << " P=2,2,2 C=8 F=simple O=0" << endl;
1255 }
1256 
1257 int main(int argc, char *argv[])
1258 {
1259  Teuchos::GlobalMPISession session(&argc, &argv);
1260  Kokkos::initialize (argc, argv);
1261  //cout << argv << endl;
1262 
1263  RCP<const Teuchos::Comm<int> > tcomm = Teuchos::DefaultComm<int>::getComm();
1264  int rank = tcomm->getRank();
1265 
1266 
1267  int numParts = -10;
1268  float imbalance = -1.03;
1269  int k = -1;
1270 
1271  string pqParts = "";
1272  int opt = 1;
1273  std::string fname = "";
1274  std::string paramFile = "";
1275 
1276 
1277  int migration_check_option = -2;
1278  int migration_all_to_all_type = -1;
1279  zscalar_t migration_imbalance_cut_off = -1.15;
1280  int migration_processor_assignment_type = -1;
1281  int migration_doMigration_type = -1;
1282  bool test_boxes = false;
1283  bool rectilinear = false;
1284 
1285  try{
1286  try {
1287  getArgVals(
1288  argc,
1289  argv,
1290  numParts,
1291  imbalance ,
1292  pqParts,
1293  opt,
1294  fname,
1295  paramFile,
1296  k,
1297  migration_check_option,
1298  migration_all_to_all_type,
1299  migration_imbalance_cut_off,
1300  migration_processor_assignment_type,
1301  migration_doMigration_type,
1302  test_boxes,
1303  rectilinear);
1304  }
1305  catch(std::string s){
1306  if(tcomm->getRank() == 0){
1307  print_usage(argv[0]);
1308  }
1309  throw s;
1310  }
1311 
1312  catch(char * s){
1313  if(tcomm->getRank() == 0){
1314  print_usage(argv[0]);
1315  }
1316  throw s;
1317  }
1318 
1319  int ierr = 0;
1320 
1321  switch (opt){
1322 
1323  case 0:
1324  ierr = testFromDataFile(tcomm,numParts, imbalance,fname,
1325  pqParts, paramFile, k,
1326  migration_check_option,
1327  migration_all_to_all_type,
1328  migration_imbalance_cut_off,
1329  migration_processor_assignment_type,
1330  migration_doMigration_type, test_boxes, rectilinear);
1331  break;
1332 #ifdef hopper_separate_test
1333  case 1:
1334  ierr = testFromSeparateDataFiles(tcomm,numParts, imbalance,fname,
1335  pqParts, paramFile, k,
1336  migration_check_option,
1337  migration_all_to_all_type,
1338  migration_imbalance_cut_off,
1339  migration_processor_assignment_type,
1340  migration_doMigration_type, test_boxes, rectilinear);
1341  break;
1342 #endif
1343  default:
1344  ierr = GeometricGenInterface(tcomm, numParts, imbalance, fname,
1345  pqParts, paramFile, k,
1346  migration_check_option,
1347  migration_all_to_all_type,
1348  migration_imbalance_cut_off,
1349  migration_processor_assignment_type,
1350  migration_doMigration_type, test_boxes, rectilinear);
1351  break;
1352  }
1353 
1354  if (rank == 0) {
1355  if (ierr == 0) std::cout << "PASS" << std::endl;
1356  else std::cout << "FAIL" << std::endl;
1357  }
1358  }
1359 
1360 
1361  catch(std::string &s){
1362  if (rank == 0)
1363  cerr << s << endl;
1364  }
1365 
1366  catch(char * s){
1367  if (rank == 0)
1368  cerr << s << endl;
1369  }
1370 
1371  Kokkos::finalize ();
1372  return 0;
1373 }
#define CATCH_EXCEPTIONS_WITH_COUNT(ierr, pp)
void print_usage(char *executable)
void readGeoGenParams(string paramFileName, Teuchos::ParameterList &geoparams, const RCP< const Teuchos::Comm< int > > &comm)
int main(int argc, char *argv[])
double zscalar_t
A simple class that can be the User template argument for an InputAdapter.
int run_pointAssign_tests(Zoltan2::PartitioningProblem< Adapter > *problem, RCP< tMVector_t > &coords)
int testFromDataFile(RCP< const Teuchos::Comm< int > > &comm, int numParts, float imbalance, std::string fname, std::string pqParts, std::string pfname, int k, int migration_check_option, int migration_all_to_all_type, zscalar_t migration_imbalance_cut_off, int migration_processor_assignment_type, int migration_doMigration_type, bool test_boxes, bool rectilinear)
static ArrayRCP< ArrayRCP< zscalar_t > > weights
int zlno_t
Defines the PartitioningSolution class.
common code used by tests
int GeometricGenInterface(RCP< const Teuchos::Comm< int > > &comm, int numParts, float imbalance, std::string paramFile, std::string pqParts, std::string pfname, int k, int migration_check_option, int migration_all_to_all_type, zscalar_t migration_imbalance_cut_off, int migration_processor_assignment_type, int migration_doMigration_type, bool test_boxes, bool rectilinear)
void print_boxAssign_result(const char *str, int dim, typename Adapter::scalar_t *lower, typename Adapter::scalar_t *upper, size_t nparts, typename Adapter::part_t *parts)
coordinateModelPartBox Class, represents the boundaries of the box which is a result of a geometric p...
Defines the XpetraMultiVectorAdapter.
const char param_comment
SparseMatrixAdapter_t::part_t part_t
void getCoords(void *data, int numGid, int numLid, int numObj, zgno_t *gids, zgno_t *lids, int dim, double *coords, int *ierr)
Defines the EvaluatePartition class.
int run_boxAssign_tests(Zoltan2::PartitioningProblem< Adapter > *problem, RCP< tMVector_t > &coords)
RCP< const Comm< int > > getComm()
Return the communicator used by the problem.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
An adapter for Xpetra::MultiVector.
#define CATCH_EXCEPTIONS_AND_RETURN(pp)
static const std::string fail
int zgno_t
string convert_to_string(char *args)
void printTimers() const
Return the communicator passed to the problem.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
string trim_right_copy(const string &s, const string &delimiters=" \\\)
PartitioningProblem sets up partitioning problems for the user.
void getArgVals(int argc, char **argv, int &numParts, float &imbalance, string &pqParts, int &opt, std::string &fname, std::string &pfname, int &k, int &migration_check_option, int &migration_all_to_all_type, zscalar_t &migration_imbalance_cut_off, int &migration_processor_assignment_type, int &migration_doMigration_type, bool &test_boxes, bool &rectilinear)
string trim_copy(const string &s, const string &delimiters=" \\\)
Tpetra::MultiVector< double, int, int > tMVector_t
int getNumWeights()
##END Predistribution functions######################//
bool getArgumentValue(string &argumentid, double &argumentValue, string argumentline)
RCP< tMVector_t > getUICoordinates()
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tMVector_t
Defines the PartitioningProblem class.
A class that computes and returns quality metrics.
string trim_left_copy(const string &s, const string &delimiters=" \\\)
Defines the BasicVectorAdapter class.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
std::string testDataFilePath(".")