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