Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ordering1.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 <fstream>
15 #include <limits>
16 #include <vector>
17 #include <Teuchos_ParameterList.hpp>
18 #include <Teuchos_RCP.hpp>
19 #include <Teuchos_CommandLineProcessor.hpp>
20 #include <Tpetra_CrsMatrix.hpp>
21 #include <Tpetra_Vector.hpp>
22 #include <MatrixMarket_Tpetra.hpp>
23 
24 using Teuchos::RCP;
25 
27 // Program to demonstrate use of Zoltan2 to order a TPetra matrix
28 // (read from a MatrixMarket file or generated by Galeri::Xpetra).
29 // Usage:
30 // a.out [--inputFile=filename] [--outputFile=outfile] [--verbose]
31 // [--x=#] [--y=#] [--z=#] [--matrix={Laplace1D,Laplace2D,Laplace3D}
32 // Karen Devine, 2011
34 
36 // Eventually want to use Teuchos unit tests to vary z2TestLO and
37 // GO. For now, we set them at compile time based on whether Tpetra
38 // is built with explicit instantiation on. (in Zoltan2_TestHelpers.hpp)
39 
40 typedef zlno_t z2TestLO;
41 typedef zgno_t z2TestGO;
43 
44 typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO> SparseMatrix;
45 typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO> Vector;
46 typedef Vector::node_type Node;
48 
49 size_t computeBandwidth(RCP<SparseMatrix> A, z2TestLO *iperm)
50 // Returns the bandwidth of the (local) permuted matrix
51 // Uses the inverse permutation calculated from the OrderingSolution
52 // if passed in, otherwise is calculating the original value.
53 {
54  z2TestLO ii, i, j, k;
55  typename SparseMatrix::local_inds_host_view_type indices;
56  typename SparseMatrix::values_host_view_type values;
57 
58  z2TestLO bw_left = 0;
59  z2TestLO bw_right = 0;
60 
61  z2TestLO n = A->getLocalNumRows();
62 
63  // Loop over rows of matrix
64  for (ii=0; ii<n; ii++) {
65  A->getLocalRowView (ii, indices, values);
66  for (k=0; k< static_cast<z2TestLO>(indices.size()); k++){
67  if (indices[k] < n){ // locally owned
68  if (iperm){
69  i = iperm[ii];
70  j = iperm[indices[k]];
71  } else {
72  i = ii;
73  j = indices[k];
74  }
75  if (j-i > bw_right)
76  bw_right = j-i;
77  if (i-j > bw_left)
78  bw_left = i-j;
79  }
80  }
81  }
82 
83  // Total bandwidth is the sum of left and right + 1
84  return (bw_left + bw_right + 1);
85 }
86 
88 int main(int narg, char** arg)
89 {
90  std::string inputFile = ""; // Matrix Market file to read
91  std::string outputFile = ""; // Output file to write
92  bool verbose = false; // Verbosity of output
93  int testReturn = 0;
94  size_t origBandwidth = 0;
95  size_t permutedBandwidth = 0;
96 
98  Tpetra::ScopeGuard tscope(&narg, &arg);
99  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
100 
101  int me = comm->getRank();
102 
103  // Read run-time options.
104  Teuchos::CommandLineProcessor cmdp (false, false);
105  cmdp.setOption("inputFile", &inputFile,
106  "Name of a Matrix Market file in the data directory; "
107  "if not specified, a matrix will be generated by Galeri.");
108  cmdp.setOption("outputFile", &outputFile,
109  "Name of file to write the permutation");
110  cmdp.setOption("verbose", "quiet", &verbose,
111  "Print messages and results.");
112  if (me == 0) std::cout << "Starting everything" << std::endl;
113 
115  // Even with cmdp option "true", I get errors for having these
116  // arguments on the command line. (On redsky build)
117  // KDDKDD Should just be warnings, right? Code should still work with these
118  // KDDKDD params in the create-a-matrix file. Better to have them where
119  // KDDKDD they are used.
120  int xdim=10;
121  int ydim=10;
122  int zdim=10;
123  std::string matrixType("Laplace3D");
124 
125  cmdp.setOption("x", &xdim,
126  "number of gridpoints in X dimension for "
127  "mesh used to generate matrix.");
128  cmdp.setOption("y", &ydim,
129  "number of gridpoints in Y dimension for "
130  "mesh used to generate matrix.");
131  cmdp.setOption("z", &zdim,
132  "number of gridpoints in Z dimension for "
133  "mesh used to generate matrix.");
134  cmdp.setOption("matrix", &matrixType,
135  "Matrix type: Laplace1D, Laplace2D, or Laplace3D");
136 
138  // Ordering options to test.
140  std::string orderMethod("rcm"); // TODO: Allow "RCM" as well
141  cmdp.setOption("order_method", &orderMethod,
142  "order_method: natural, random, rcm, sorted_degree");
143 
144  std::string orderMethodType("local"); // TODO: Allow "LOCAL" as well
145  cmdp.setOption("order_method_type", &orderMethodType,
146  "local or global or both");
147 
149  cmdp.parse(narg, arg);
150 
151 
152  RCP<UserInputForTests> uinput;
153 
154  if (inputFile != ""){ // Input file specified; read a matrix
155  uinput = rcp(new UserInputForTests(testDataFilePath, inputFile, comm, true));
156  }
157  else // Let Galeri generate a matrix
158 
159  uinput = rcp(new UserInputForTests(xdim, ydim, zdim, matrixType, comm, true, true));
160 
161  RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
162 
163  if (me == 0)
164  std::cout << "NumRows = " << origMatrix->getGlobalNumRows() << std::endl
165  << "NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
166  << "NumProcs = " << comm->getSize() << std::endl;
167 
169  // Currently Not Used
170  /*
171  RCP<Vector> origVector, origProd;
172  origProd = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
173  origMatrix->getRangeMap());
174  origVector = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
175  origMatrix->getDomainMap());
176  origVector->randomize();
177  */
178 
180  Teuchos::ParameterList params;
181  params.set("order_method", orderMethod);
182  params.set("order_method_type", orderMethodType);
183 
185  SparseMatrixAdapter adapter(origMatrix);
186 
188  try
189  {
190  Zoltan2::OrderingProblem<SparseMatrixAdapter> problem(&adapter, &params);
191  if (me == 0) std::cout << "Going to solve" << std::endl;
192  problem.solve();
193 
195 
197  problem.getLocalOrderingSolution();
198 
199  if (me == 0) std::cout << "Going to get results" << std::endl;
200  // Check that the solution is really a permutation
201 
202  z2TestLO * perm = soln->getPermutationView();
203 
204  if (outputFile != "") {
205  std::ofstream permFile;
206 
207  // Write permutation (0-based) to file
208  // each process writes local perm to a separate file
209  //std::string fname = outputFile + "." + me;
210  std::stringstream fname;
211  fname << outputFile << "." << comm->getSize() << "." << me;
212  permFile.open(fname.str().c_str());
213  size_t checkLength = soln->getPermutationSize();
214  for (size_t i=0; i<checkLength; i++){
215  permFile << " " << perm[i] << std::endl;
216  }
217  permFile.close();
218 
219  }
220 
221  if (me == 0) std::cout << "Verifying results " << std::endl;
222 
223  // Verify that checkPerm is a permutation
224  testReturn = soln->validatePerm();
225 
226  // Compute original bandwidth
227  origBandwidth = computeBandwidth(origMatrix, nullptr);
228 
229  // Compute permuted bandwidth
230  z2TestLO * iperm = soln->getPermutationView(true);
231  permutedBandwidth = computeBandwidth(origMatrix, iperm);
232 
233  } catch (std::exception &e){
234  std::cout << "Exception caught in ordering" << std::endl;
235  std::cout << e.what() << std::endl;
236  std::cout << "FAIL" << std::endl;
237  return 0;
238  }
239 
240  if (testReturn)
241  std::cout << me << ": Not a valid permutation" << std::endl;
242 
243  int gTestReturn;
244  Teuchos::reduceAll<int,int>(*comm, Teuchos::EReductionType::REDUCE_MAX, 1,
245  &testReturn, &gTestReturn);
246 
247  int increasedBandwidth = (permutedBandwidth > origBandwidth);
248  if (increasedBandwidth)
249  std::cout << me << ": Bandwidth increased: original "
250  << origBandwidth << " < " << permutedBandwidth << " permuted "
251  << std::endl;
252  else
253  std::cout << me << ": Bandwidth not increased: original "
254  << origBandwidth << " >= " << permutedBandwidth << " permuted "
255  << std::endl;
256 
257  int gIncreasedBandwidth;
258  Teuchos::reduceAll<int,int>(*comm, Teuchos::EReductionType::REDUCE_MAX, 1,
259  &increasedBandwidth, &gIncreasedBandwidth);
260 
261  if (me == 0) {
262  if (gTestReturn)
263  std::cout << "Solution is not a permutation; FAIL" << std::endl;
264  else if (gIncreasedBandwidth && (orderMethod == "rcm"))
265  std::cout << "Bandwidth was increased; FAIL" << std::endl;
266  else
267  std::cout << "PASS" << std::endl;
268 
269  }
270 
271 }
272 
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
size_t computeBandwidth(RCP< SparseMatrix > A, z2TestLO *iperm)
Definition: ordering1.cpp:49
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
Definition: coloring1.cpp:45
list fname
Begin.
Definition: validXML.py:19
OrderingProblem sets up ordering problems for the user.
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.
int validatePerm()
returns 0 if permutation is valid, negative if invalid.
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(".")