Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
orderingMetis.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Zoltan2: A package of combinatorial algorithms for scientific computing
4 //
5 // Copyright 2012 NTESS and the Zoltan2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
12 #include <Zoltan2_TestHelpers.hpp>
13 #include <iostream>
14 #include <limits>
15 #include <vector>
16 #include <Teuchos_ParameterList.hpp>
17 #include <Teuchos_RCP.hpp>
18 #include <Teuchos_CommandLineProcessor.hpp>
19 #include <Tpetra_CrsMatrix.hpp>
20 #include <Tpetra_Vector.hpp>
21 #include <MatrixMarket_Tpetra.hpp>
22 
23 //#include <Zoltan2_Memory.hpp> KDD User app wouldn't include our memory mgr.
24 
25 using Teuchos::RCP;
26 
28 // Program to demonstrate use of Zoltan2 to order a TPetra matrix
29 // (read from a MatrixMarket file or generated by Galeri::Xpetra).
30 // Usage:
31 // a.out [--inputFile=filename] [--outputFile=outfile] [--verbose]
32 // [--x=#] [--y=#] [--z=#] [--matrix={Laplace1D,Laplace2D,Laplace3D}
33 // Karen Devine, 2011
35 
37 // Eventually want to use Teuchos unit tests to vary z2TestLO and
38 // GO. For now, we set them at compile time based on whether Tpetra
39 // is built with explicit instantiation on. (in Zoltan2_TestHelpers.hpp)
40 
41 typedef zlno_t z2TestLO;
42 typedef zgno_t z2TestGO;
44 
45 typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO> SparseMatrix;
46 typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO> Vector;
47 typedef Vector::node_type Node;
48 
50 
51 #define epsilon 0.00000001
52 
53 int validatePerm(size_t n, z2TestLO *perm)
54 // returns 0 if permutation is valid
55 {
56  std::vector<int> count(n);
57  int status = 0;
58 
59  for (size_t i=0; i<n; i++)
60  count[i]=0;
61 
62  for (size_t i=0; i<n; i++){
63  if (perm[i] < 0 || static_cast<size_t> (perm[i]) >= n)
64  status = -1;
65  else
66  count[perm[i]]++;
67  }
68 
69  // Each index should occur exactly once (count==1)
70  for (size_t i=0; i<n; i++){
71  if (count[i] != 1){
72  status = -2;
73  break;
74  }
75  }
76 
77  return status;
78 }
79 
81 int main(int narg, char** arg)
82 {
83  std::string inputFile = ""; // Matrix Market file to read
84  std::string outputFile = ""; // Matrix Market file to write
85  bool verbose = false; // Verbosity of output
86  int testReturn = 0;
87 
88  // Initialize Tpetra and get default communicator.
89  Tpetra::ScopeGuard tscope(&narg, &arg);
90  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
91  int me = comm->getRank();
92 
93  // Read run-time options.
94  Teuchos::CommandLineProcessor cmdp (false, false);
95  cmdp.setOption("inputFile", &inputFile,
96  "Name of a Matrix Market file in the data directory; "
97  "if not specified, a matrix will be generated by MueLu.");
98  cmdp.setOption("outputFile", &outputFile,
99  "Name of the Matrix Market sparse matrix file to write, "
100  "echoing the input/generated matrix.");
101  cmdp.setOption("verbose", "quiet", &verbose,
102  "Print messages and results.");
103 
105  // Even with cmdp option "true", I get errors for having these
106  // arguments on the command line. (On redsky build)
107  // KDDKDD Should just be warnings, right? Code should still work with these
108  // KDDKDD params in the create-a-matrix file. Better to have them where
109  // KDDKDD they are used.
110  int xdim=10;
111  int ydim=10;
112  int zdim=10;
113  std::string matrixType("Laplace3D");
114 
115  cmdp.setOption("x", &xdim,
116  "number of gridpoints in X dimension for "
117  "mesh used to generate matrix.");
118  cmdp.setOption("y", &ydim,
119  "number of gridpoints in Y dimension for "
120  "mesh used to generate matrix.");
121  cmdp.setOption("z", &zdim,
122  "number of gridpoints in Z dimension for "
123  "mesh used to generate matrix.");
124  cmdp.setOption("matrix", &matrixType,
125  "Matrix type: Laplace1D, Laplace2D, or Laplace3D");
127 
128  cmdp.parse(narg, arg);
129 
130 
131  RCP<UserInputForTests> uinput;
132 
133  if (inputFile != "") // Input file specified; read a matrix
134  uinput = rcp(new UserInputForTests(testDataFilePath, inputFile, comm, true));
135 
136  else // Let MueLu generate a matrix
137  uinput = rcp(new UserInputForTests(xdim, ydim, zdim, matrixType, comm, true, true));
138 
139  RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
140 
141  if (outputFile != "") {
142  // Just a sanity check.
143  Tpetra::MatrixMarket::Writer<SparseMatrix>::writeSparseFile(outputFile,
144  origMatrix, verbose);
145  }
146 
147  if (me == 0)
148  std::cout << "NumRows = " << origMatrix->getGlobalNumRows() << std::endl
149  << "NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
150  << "NumProcs = " << comm->getSize() << std::endl;
151 
153  RCP<Vector> origVector, origProd;
154  origProd = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
155  origMatrix->getRangeMap());
156  origVector = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
157  origMatrix->getDomainMap());
158  origVector->randomize();
159 
161  Teuchos::ParameterList params;
163  size_t checkLength;
164  z2TestLO *checkPerm;
165 
167  SparseMatrixAdapter adapter(origMatrix);
168 
169  params.set("order_method", "metis");
170  params.set("order_method_type", "local");
171 
173  try
174  {
175  Zoltan2::OrderingProblem<SparseMatrixAdapter> problem(&adapter, &params);
176  problem.solve();
177 
179  problem.getLocalOrderingSolution();
180 
181  // Check that the solution is really a permutation
182  checkLength = soln->getPermutationSize();
183  checkPerm = soln->getPermutationView();
184 
185  for (size_t ii = 0; ii < checkLength; ii++)
186  std::cout << checkPerm[ii] << " ";
187  std::cout << std::endl;
188  // Verify that checkPerm is a permutation
189  testReturn = validatePerm(checkLength, checkPerm);
190 
191  } catch (std::exception &e){
192 #ifdef HAVE_ZOLTAN2_METIS
193  // AMD is defined and still got an exception.
194  std::cout << "Exception from METIS Algorithm" << std::endl;
195  std::cout << "FAIL" << std::endl;
196  return 0;
197 #else
198  std::cout << "METIS is not enabled. METIS Algorithm threw an exception."
199  << std::endl;
200  std::cout << "PASS" << std::endl;
201  return 0;
202 #endif
203  }
204 
205  if (me == 0) {
206  if (testReturn)
207  std::cout << "Solution is not a permutation; FAIL" << std::endl;
208  else
209  std::cout << "PASS" << std::endl;
210  }
211  return 0;
212 }
213 
zgno_t z2TestGO
Definition: coloring1.cpp:41
zlno_t z2TestLO
Definition: coloring1.cpp:40
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
Tpetra::CrsMatrix< z2TestScalar, z2TestLO, z2TestGO > SparseMatrix
Definition: coloring1.cpp:44
int main(int narg, char **arg)
Definition: coloring1.cpp:164
common code used by tests
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
Definition: coloring1.cpp:45
OrderingProblem sets up ordering problems for the user.
int validatePerm(size_t n, z2TestLO *perm)
Definition: orderingAMD.cpp:53
Defines the XpetraCrsMatrixAdapter class.
size_t getPermutationSize() const
Get (local) size of permutation.
lno_t * getPermutationView(bool inverse=false) const
Get pointer to permutation. If inverse = true, return inverse permutation. By default, perm[i] is where new index i can be found in the old ordering. When inverse==true, perm[i] is where old index i can be found in the new ordering.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Tpetra::Map::local_ordinal_type zlno_t
Defines the OrderingProblem class.
Zoltan2::XpetraCrsMatrixAdapter< SparseMatrix > SparseMatrixAdapter
Definition: coloring1.cpp:50
zscalar_t z2TestScalar
Definition: coloring1.cpp:42
float zscalar_t
Tpetra::Map::global_ordinal_type zgno_t
Vector::node_type Node
Definition: coloring1.cpp:46
LocalOrderingSolution< lno_t > * getLocalOrderingSolution()
Get the local ordering solution to the problem.
std::string testDataFilePath(".")