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>
44 typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO>
SparseMatrix;
45 typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO>
Vector;
46 typedef Vector::node_type
Node;
55 typename SparseMatrix::local_inds_host_view_type indices;
56 typename SparseMatrix::values_host_view_type values;
64 for (ii=0; ii<n; ii++) {
65 A->getLocalRowView (ii, indices, values);
66 for (k=0; k< static_cast<z2TestLO>(indices.size()); k++){
70 j = iperm[indices[k]];
84 return (bw_left + bw_right + 1);
88 int main(
int narg,
char** arg)
90 std::string inputFile =
"";
91 std::string outputFile =
"";
94 size_t origBandwidth = 0;
95 size_t permutedBandwidth = 0;
98 Tpetra::ScopeGuard tscope(&narg, &arg);
99 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
101 int me = comm->getRank();
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;
123 std::string matrixType(
"Laplace3D");
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");
140 std::string orderMethod(
"rcm");
141 cmdp.setOption(
"order_method", &orderMethod,
142 "order_method: natural, random, rcm, sorted_degree");
144 std::string orderMethodType(
"local");
145 cmdp.setOption(
"order_method_type", &orderMethodType,
146 "local or global or both");
149 cmdp.parse(narg, arg);
152 RCP<UserInputForTests> uinput;
154 if (inputFile !=
""){
159 uinput = rcp(
new UserInputForTests(xdim, ydim, zdim, matrixType, comm,
true,
true));
161 RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
164 std::cout <<
"NumRows = " << origMatrix->getGlobalNumRows() << std::endl
165 <<
"NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
166 <<
"NumProcs = " << comm->getSize() << std::endl;
180 Teuchos::ParameterList params;
181 params.set(
"order_method", orderMethod);
182 params.set(
"order_method_type", orderMethodType);
191 if (me == 0) std::cout <<
"Going to solve" << std::endl;
199 if (me == 0) std::cout <<
"Going to get results" << std::endl;
204 if (outputFile !=
"") {
205 std::ofstream permFile;
210 std::stringstream
fname;
211 fname << outputFile <<
"." << comm->getSize() <<
"." << me;
212 permFile.open(fname.str().c_str());
214 for (
size_t i=0; i<checkLength; i++){
215 permFile <<
" " << perm[i] << std::endl;
221 if (me == 0) std::cout <<
"Verifying results " << std::endl;
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;
241 std::cout << me <<
": Not a valid permutation" << std::endl;
244 Teuchos::reduceAll<int,int>(*comm, Teuchos::EReductionType::REDUCE_MAX, 1,
245 &testReturn, &gTestReturn);
247 int increasedBandwidth = (permutedBandwidth > origBandwidth);
248 if (increasedBandwidth)
249 std::cout << me <<
": Bandwidth increased: original "
250 << origBandwidth <<
" < " << permutedBandwidth <<
" permuted "
253 std::cout << me <<
": Bandwidth not increased: original "
254 << origBandwidth <<
" >= " << permutedBandwidth <<
" permuted "
257 int gIncreasedBandwidth;
258 Teuchos::reduceAll<int,int>(*comm, Teuchos::EReductionType::REDUCE_MAX, 1,
259 &increasedBandwidth, &gIncreasedBandwidth);
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;
267 std::cout <<
"PASS" << std::endl;
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
Tpetra::CrsMatrix< z2TestScalar, z2TestLO, z2TestGO > SparseMatrix
int main(int narg, char **arg)
common code used by tests
size_t computeBandwidth(RCP< SparseMatrix > A, z2TestLO *iperm)
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
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
Tpetra::Map::global_ordinal_type zgno_t
LocalOrderingSolution< lno_t > * getLocalOrderingSolution()
Get the local ordering solution to the problem.
std::string testDataFilePath(".")