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>
45 typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO>
SparseMatrix;
46 typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO>
Vector;
47 typedef Vector::node_type
Node;
51 #define epsilon 0.00000001
56 std::vector<int> count(n);
59 for (
size_t i=0; i<n; i++)
62 for (
size_t i=0; i<n; i++){
63 if (perm[i] < 0 || static_cast<size_t> (perm[i]) >= n)
70 for (
size_t i=0; i<n; i++){
81 int main(
int narg,
char** arg)
83 std::string inputFile =
"";
84 std::string outputFile =
"";
89 Tpetra::ScopeGuard tscope(&narg, &arg);
90 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
91 int me = comm->getRank();
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.");
113 std::string matrixType(
"Laplace3D");
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");
128 cmdp.parse(narg, arg);
131 RCP<UserInputForTests> uinput;
137 uinput = rcp(
new UserInputForTests(xdim, ydim, zdim, matrixType, comm,
true,
true));
139 RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
141 if (outputFile !=
"") {
143 Tpetra::MatrixMarket::Writer<SparseMatrix>::writeSparseFile(outputFile,
144 origMatrix, verbose);
148 std::cout <<
"NumRows = " << origMatrix->getGlobalNumRows() << std::endl
149 <<
"NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
150 <<
"NumProcs = " << comm->getSize() << std::endl;
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();
161 Teuchos::ParameterList params;
169 params.set(
"order_method",
"metis");
170 params.set(
"order_method_type",
"local");
185 for (
size_t ii = 0; ii < checkLength; ii++)
186 std::cout << checkPerm[ii] <<
" ";
187 std::cout << std::endl;
191 }
catch (std::exception &e){
192 #ifdef HAVE_ZOLTAN2_METIS
194 std::cout <<
"Exception from METIS Algorithm" << std::endl;
195 std::cout <<
"FAIL" << std::endl;
198 std::cout <<
"METIS is not enabled. METIS Algorithm threw an exception."
200 std::cout <<
"PASS" << std::endl;
207 std::cout <<
"Solution is not a permutation; FAIL" << std::endl;
209 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
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
OrderingProblem sets up ordering problems for the user.
int validatePerm(size_t n, z2TestLO *perm)
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
Tpetra::Map::global_ordinal_type zgno_t
LocalOrderingSolution< lno_t > * getLocalOrderingSolution()
Get the local ordering solution to the problem.
std::string testDataFilePath(".")