Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
partition_sarma.cpp
Go to the documentation of this file.
1 /*
2  * Created by mbenlioglu on Nov 10, 2020.
3  */
4 
5 #include <Zoltan2_config.h>
6 
11 
12 using Teuchos::RCP;
13 using Teuchos::rcp;
14 
16 
17 #ifdef HAVE_ZOLTAN2_SARMA
18 
19 struct Parameters {
20  std::string alg = "pal";
21  sarma::Order order_type = sarma::Order::NAT;
22  std::string order_str = "nat";
23  Ordinal row_parts = 8, col_parts = 0;
24  Value max_load = 0;
25  int seed = 2147483647;
26  double sparsify = 1.0;
27  bool triangular = false, use_data = false;
28 };
29 
30 const static std::vector<std::pair<Parameters, std::vector<Ordinal> > > expected = {
31  {{.alg="opal"}, {0, 48, 94, 143, 175, 218, 257, 295, 332, 0, 48, 94, 143, 175, 218, 257, 295, 332}}
32 };
33 
34 auto testFromFile(const RCP<const Teuchos::Comm<int> > &comm, int nparts, std::string &filename, bool doRemap,
35  const Parameters sarma) {
36  int me = comm->getRank();
37 
38  UserInputForTests userInputForTests(testDataFilePath, filename, comm, false, true);
39  RCP<tcrsMatrix_t> matrix = userInputForTests.getUITpetraCrsMatrix();
40  RCP<const tcrsMatrix_t> matrixConst = rcp_const_cast<const tcrsMatrix_t>(matrix);
41 
42 
43  xCM_tCM_t matrixAdapter(matrixConst, 0);
44 
45  #ifdef HAVE_ZOLTAN2_MPI
46  double begin = MPI_Wtime();
47  #endif
48 
49  Teuchos::ParameterList params("test params");
50  params.set("debug_level", "basic_status");
51  params.set("num_global_parts", nparts);
52  params.set("algorithm", "sarma");
53  if (doRemap) params.set("remap_parts", true);
54  Teuchos::ParameterList &zparams = params.sublist("zoltan_parameters", false);
55  zparams.set("DEBUG_LEVEL", "0");
56 
57  // Sarma params
58  Teuchos::ParameterList &sparams = params.sublist("sarma_parameters", false);
59  sparams.set("alg", sarma.alg);
60  sparams.set("order", sarma.order_str);
61  sparams.set("row_parts", sarma.row_parts);
62  sparams.set("col_parts", sarma.col_parts);
63  sparams.set("max_load", sarma.max_load);
64  sparams.set("sparsify", sarma.sparsify);
65  sparams.set("triangular", sarma.triangular);
66  sparams.set("use_data", sarma.use_data);
67  sparams.set("seed", sarma.seed);
68 
69  #ifdef HAVE_ZOLTAN2_MPI
70  Zoltan2::PartitioningProblem<xCM_tCM_t> problem(&matrixAdapter, &params, comm);
71  #else
72  Zoltan2::PartitioningProblem<xCM_tCM_t> problem(&matrixAdapter, &params);
73  #endif
74  if (me == 0) std::cout << "Problem constructed" << std::endl;
75 
76  problem.solve();
77 
78  #ifdef HAVE_ZOLTAN2_MPI
79  if (me == 0)
80  std::cout << "Run-time:\t" << MPI_Wtime() - begin << std::endl;
81  #endif
82  return problem.getSolution();
83 }
84 #endif
85 
86 int main(int argc, char *argv[]) {
87  #ifdef HAVE_ZOLTAN2_SARMA
88  Tpetra::ScopeGuard tscope(&argc, &argv);
89  RCP<const Teuchos::Comm<int> > tcomm = Tpetra::getDefaultComm();
90 
91  Parameters params;
92 
93  int nParts = 8;
94  bool doRemap = false;
95  std::string filename = "USAir97";
96 
97  // Run-time options
98 
99  for (auto test : expected){
100 
101  auto zsoln = testFromFile(tcomm, nParts, filename, doRemap, test.first);
102 
103  // compare zoltan vs raw
104 
105  const int *zparts = zsoln.getPartListView();
106  for (unsigned i = 0; i < test.second.size(); ++i) {
107  if (zparts[i] != (int) test.second[i]) {
108  std::cout << "FAIL" << std::endl;
109  return EXIT_FAILURE;
110  }
111  }
112  }
113  #endif
114  std::cout << "PASS" << std::endl;
115 
116  return EXIT_SUCCESS;
117 }
#define nParts
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
int main(int narg, char **arg)
Definition: coloring1.cpp:199
Defines the PartitioningSolution class.
common code used by tests
Defines the XpetraCrsMatrixAdapter class.
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > myTypes_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tcrsMatrix_t
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
std::string testDataFilePath(".")