44 typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO>
SparseMatrix;
46 typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO>
Vector;
47 typedef Vector::node_type
Node;
57 typedef Tpetra::Vector<int, z2TestLO, z2TestGO>
IntVector;
60 #define epsilon 0.00000001
64 int main(
int narg,
char** arg)
66 std::string inputFile =
"";
72 std::string method =
"scotch";
74 bool distributeInput =
true;
81 Tpetra::ScopeGuard tscope(&narg, &arg);
82 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
83 int me = comm->getRank();
86 Teuchos::CommandLineProcessor cmdp (
false,
false);
87 cmdp.setOption(
"inputPath", &inputPath,
88 "Path to the MatrixMarket or Zoltan file to be read; ",
true);
90 cmdp.setOption(
"inputFile", &inputFile,
91 "Name of the Matrix Market or Zoltan file to read; ",
true);
97 cmdp.setOption(
"verbose",
"quiet", &verbose,
98 "Print messages and results.");
104 cmdp.parse(narg, arg);
106 RCP<UserInputForTests> uinput;
109 uinput = rcp(
new UserInputForTests(inputPath, inputFile, comm,
true, distributeInput));
111 RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
113 if (origMatrix->getGlobalNumRows() < 40)
115 Teuchos::FancyOStream out(Teuchos::rcp(&std::cout,
false));
116 origMatrix->describe(out, Teuchos::VERB_EXTREME);
122 std::cout <<
"NumRows = " << origMatrix->getGlobalNumRows() << std::endl
123 <<
"NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
124 <<
"NumProcs = " << comm->getSize() << std::endl
125 <<
"NumLocalRows (rank 0) = " << origMatrix->getLocalNumRows() << std::endl;
139 Teuchos::ParameterList params;
141 params.set(
"partitioning_approach",
"partition");
161 if (me == 0) std::cout <<
"Calling solve() " << std::endl;
165 if (me == 0) std::cout <<
"Done solve() " << std::endl;
168 catch (std::runtime_error &e)
170 std::cout <<
"Runtime exception returned from solve(): " << e.what();
171 if (!strncmp(e.what(),
"BUILD ERROR", 11)) {
173 std::cout <<
" PASS" << std::endl;
178 std::cout <<
" FAIL" << std::endl;
182 catch (std::logic_error &e) {
183 std::cout <<
"Logic exception returned from solve(): " << e.what()
184 <<
" FAIL" << std::endl;
187 catch (std::bad_alloc &e) {
188 std::cout <<
"Bad_alloc exception returned from solve(): " << e.what()
189 <<
" FAIL" << std::endl;
192 catch (std::exception &e) {
193 std::cout <<
"Unknown exception returned from solve(). " << e.what()
194 <<
" FAIL" << std::endl;
281 std::cout <<
"Finished" << std::endl;
MatrixPartitioningProblem sets up partitioning problems for the user.
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
Defines the XpetraCrsMatrixAdapter class.
Defines the MatrixPartitioningProblem class.
An adapter for Xpetra::MultiVector.
Tpetra::Map::local_ordinal_type zlno_t
Zoltan2::XpetraMultiVectorAdapter< IntVector > IntVectorAdapter
Tpetra::CrsGraph< z2TestLO, z2TestGO > SparseGraph
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Zoltan2::XpetraCrsMatrixAdapter< SparseMatrix > SparseMatrixAdapter
Tpetra::Vector< int, z2TestLO, z2TestGO > IntVector
Zoltan2::XpetraMultiVectorAdapter< Vector > MultiVectorAdapter
Tpetra::Map::global_ordinal_type zgno_t
std::string testDataFilePath(".")