52 #include <Teuchos_ParameterList.hpp>
53 #include <Teuchos_RCP.hpp>
54 #include <Teuchos_FancyOStream.hpp>
55 #include <Teuchos_CommandLineProcessor.hpp>
56 #include <Tpetra_CrsMatrix.hpp>
57 #include <Tpetra_Vector.hpp>
58 #include <MatrixMarket_Tpetra.hpp>
80 typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO>
SparseMatrix;
82 typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO>
Vector;
83 typedef Vector::node_type
Node;
91 typedef Tpetra::Vector<int, z2TestLO, z2TestGO>
IntVector;
94 #define epsilon 0.00000001
98 int main(
int narg,
char** arg)
100 std::string inputFile =
"";
101 std::string outputFile =
"";
103 std::string method =
"scotch";
104 bool verbose =
false;
105 bool distributeInput =
true;
106 bool haveFailure =
false;
112 Tpetra::ScopeGuard tscope(&narg, &arg);
113 RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
114 int me = comm->getRank();
117 Teuchos::CommandLineProcessor cmdp (
false,
false);
118 cmdp.setOption(
"inputPath", &inputPath,
119 "Path to the MatrixMarket or Zoltan file to be read; "
120 "if not specified, a default path will be used.");
121 cmdp.setOption(
"inputFile", &inputFile,
122 "Name of the Matrix Market or Zoltan file to read; "
123 "if not specified, a matrix will be generated by MueLu.");
124 cmdp.setOption(
"outputFile", &outputFile,
125 "Name of the Matrix Market sparse matrix file to write, "
126 "echoing the input/generated matrix.");
127 cmdp.setOption(
"method", &method,
128 "Partitioning method to use: scotch or parmetis.");
129 cmdp.setOption(
"vertexWeights", &nVwgts,
130 "Number of weights to generate for each vertex");
131 cmdp.setOption(
"edgeWeights", &nEwgts,
132 "Number of weights to generate for each edge");
133 cmdp.setOption(
"verbose",
"quiet", &verbose,
134 "Print messages and results.");
135 cmdp.setOption(
"distribute",
"no-distribute", &distributeInput,
136 "indicate whether or not to distribute "
137 "input across the communicator");
148 std::string matrixType(
"Laplace3D");
150 cmdp.setOption(
"x", &xdim,
151 "number of gridpoints in X dimension for "
152 "mesh used to generate matrix.");
153 cmdp.setOption(
"y", &ydim,
154 "number of gridpoints in Y dimension for "
155 "mesh used to generate matrix.");
156 cmdp.setOption(
"z", &zdim,
157 "number of gridpoints in Z dimension for "
158 "mesh used to generate matrix.");
159 cmdp.setOption(
"matrix", &matrixType,
160 "Matrix type: Laplace1D, Laplace2D, or Laplace3D");
163 cmdp.parse(narg, arg);
165 RCP<UserInputForTests> uinput;
169 true, distributeInput));
173 true, distributeInput));
175 RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
177 if (origMatrix->getGlobalNumRows() < 40) {
178 Teuchos::FancyOStream out(Teuchos::rcp(&std::cout,
false));
179 origMatrix->describe(out, Teuchos::VERB_EXTREME);
183 if (outputFile !=
"") {
185 Tpetra::MatrixMarket::Writer<SparseMatrix>::writeSparseFile(outputFile,
186 origMatrix, verbose);
190 std::cout <<
"NumRows = " << origMatrix->getGlobalNumRows() << std::endl
191 <<
"NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
192 <<
"NumProcs = " << comm->getSize() << std::endl
193 <<
"NumLocalRows (rank 0) = " << origMatrix->getNodeNumRows() << std::endl;
196 RCP<Vector> origVector, origProd;
197 origProd = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
198 origMatrix->getRangeMap());
199 origVector = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
200 origMatrix->getDomainMap());
201 origVector->randomize();
204 Teuchos::ParameterList params;
206 params.set(
"partitioning_approach",
"partition");
207 params.set(
"algorithm", method);
219 size_t nrows = origMatrix->getNodeNumRows();
222 for (
size_t i = 0; i < nrows; i++) {
223 size_t idx = i * nVwgts;
224 vwgts[idx] =
zscalar_t(origMatrix->getRowMap()->getGlobalElement(i))
226 for (
int j = 1; j < nVwgts; j++) vwgts[idx+j] = 1.;
228 for (
int j = 0; j < nVwgts; j++) {
229 if (j !=
NNZ_IDX) adapter.setVertexWeights(&vwgts[j], nVwgts, j);
230 else adapter.setVertexWeightIsDegree(
NNZ_IDX);
237 size_t nnz = origMatrix->getNodeNumEntries();
239 size_t nrows = origMatrix->getNodeNumRows();
240 size_t maxnzrow = origMatrix->getNodeMaxNumRowEntries();
243 Array<z2TestGO> egids(maxnzrow);
244 Array<zscalar_t> evals(maxnzrow);
245 for (
size_t i = 0; i < nrows; i++) {
247 z2TestGO gid = origMatrix->getRowMap()->getGlobalElement(i);
248 origMatrix->getGlobalRowCopy(gid, egids(), evals(), nnzinrow);
249 for (
size_t k = 0; k < nnzinrow; k++) {
250 ewgts[cnt] = (gid < egids[k] ? gid : egids[k]);
251 if (nEwgts > 1) ewgts[cnt+nnz] = (gid < egids[k] ? egids[k] : gid);
252 for (
int j = 2; j < nEwgts; j++) ewgts[cnt+nnz*j] = 1.;
256 for (
int j = 0; j < nEwgts; j++) {
257 adapter.setEdgeWeights(&ewgts[j*nnz], 1, j);
267 if (me == 0) std::cout <<
"Calling solve() " << std::endl;
271 if (me == 0) std::cout <<
"Done solve() " << std::endl;
273 catch (std::runtime_error &e) {
276 std::cout <<
"Runtime exception returned from solve(): " << e.what();
277 if (!strncmp(e.what(),
"BUILD ERROR", 11)) {
279 std::cout <<
" PASS" << std::endl;
284 std::cout <<
" FAIL" << std::endl;
288 catch (std::logic_error &e) {
291 std::cout <<
"Logic exception returned from solve(): " << e.what()
292 <<
" FAIL" << std::endl;
295 catch (std::bad_alloc &e) {
298 std::cout <<
"Bad_alloc exception returned from solve(): " << e.what()
299 <<
" FAIL" << std::endl;
302 catch (std::exception &e) {
305 std::cout <<
"Unknown exception returned from solve(). " << e.what()
306 <<
" FAIL" << std::endl;
312 size_t checkNparts = comm->getSize();
314 size_t checkLength = origMatrix->getNodeNumRows();
318 size_t *countPerPart =
new size_t[checkNparts];
319 size_t *globalCountPerPart =
new size_t[checkNparts];
322 for (
size_t i = 0; i < checkNparts; i++) countPerPart[i] = 0;
323 for (
size_t i = 0; i < checkNparts * nVwgts; i++) wtPerPart[i] = 0.;
325 for (
size_t i = 0; i < checkLength; i++) {
326 if (
size_t(checkParts[i]) >= checkNparts)
327 std::cout <<
"Invalid Part " << checkParts[i] <<
": FAIL" << std::endl;
328 countPerPart[checkParts[i]]++;
329 for (
int j = 0; j < nVwgts; j++) {
331 wtPerPart[checkParts[i]*nVwgts+j] += vwgts[i*nVwgts+j];
333 wtPerPart[checkParts[i]*nVwgts+j] += origMatrix->getNumEntriesInLocalRow(i);
336 Teuchos::reduceAll<int, size_t>(*comm, Teuchos::REDUCE_SUM, checkNparts,
337 countPerPart, globalCountPerPart);
338 Teuchos::reduceAll<int, zscalar_t>(*comm, Teuchos::REDUCE_SUM,
340 wtPerPart, globalWtPerPart);
342 size_t min = std::numeric_limits<std::size_t>::max();
345 size_t minrank = 0, maxrank = 0;
346 for (
size_t i = 0; i < checkNparts; i++) {
347 if (globalCountPerPart[i] < min) {min = globalCountPerPart[i]; minrank = i;}
348 if (globalCountPerPart[i] > max) {max = globalCountPerPart[i]; maxrank = i;}
349 sum += globalCountPerPart[i];
353 float avg = (float) sum / (
float) checkNparts;
354 std::cout <<
"Minimum count: " << min <<
" on rank " << minrank << std::endl;
355 std::cout <<
"Maximum count: " << max <<
" on rank " << maxrank << std::endl;
356 std::cout <<
"Average count: " << avg << std::endl;
357 std::cout <<
"Total count: " << sum
358 << (sum != origMatrix->getGlobalNumRows()
359 ?
"Work was lost; FAIL"
362 std::cout <<
"Imbalance: " << max / avg << std::endl;
364 std::vector<zscalar_t> minwt(nVwgts, std::numeric_limits<zscalar_t>::max());
365 std::vector<zscalar_t> maxwt(nVwgts, 0.);
366 std::vector<zscalar_t> sumwt(nVwgts, 0.);
367 for (
size_t i = 0; i < checkNparts; i++) {
368 for (
int j = 0; j < nVwgts; j++) {
369 size_t idx = i*nVwgts+j;
370 if (globalWtPerPart[idx] < minwt[j]) minwt[j] = globalWtPerPart[idx];
371 if (globalWtPerPart[idx] > maxwt[j]) maxwt[j] = globalWtPerPart[idx];
372 sumwt[j] += globalWtPerPart[idx];
375 for (
int j = 0; j < nVwgts; j++) {
376 float avgwt = (float) sumwt[j] / (
float) checkNparts;
377 std::cout << std::endl;
378 std::cout <<
"Minimum weight[" << j <<
"]: " << minwt[j] << std::endl;
379 std::cout <<
"Maximum weight[" << j <<
"]: " << maxwt[j] << std::endl;
380 std::cout <<
"Average weight[" << j <<
"]: " << avgwt << std::endl;
381 std::cout <<
"Imbalance: " << maxwt[j] / avgwt << std::endl;
386 delete [] countPerPart;
388 delete [] globalCountPerPart;
389 delete [] globalWtPerPart;
395 if (me == 0) std::cout <<
"Redistributing matrix..." << std::endl;
400 if (redistribMatrix->getGlobalNumRows() < 40) {
401 Teuchos::FancyOStream out(Teuchos::rcp(&std::cout,
false));
402 redistribMatrix->describe(out, Teuchos::VERB_EXTREME);
405 if (me == 0) std::cout <<
"Redistributing vectors..." << std::endl;
413 RCP<Vector> redistribProd;
414 redistribProd = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
415 redistribMatrix->getRangeMap());
419 RCP<IntVector> origIntVec;
421 origIntVec = Tpetra::createVector<int,z2TestLO,z2TestGO>(
422 origMatrix->getRangeMap());
423 for (
size_t i = 0; i < origIntVec->getLocalLength(); i++)
424 origIntVec->replaceLocalValue(i, me);
429 int origIntNorm = origIntVec->norm1();
430 int redistIntNorm = redistIntVec->norm1();
431 if (me == 0) std::cout <<
"IntegerVectorTest: " << origIntNorm <<
" == "
432 << redistIntNorm <<
" ?";
433 if (origIntNorm != redistIntNorm) {
434 if (me == 0) std::cout <<
" FAIL" << std::endl;
437 else if (me == 0) std::cout <<
" OK" << std::endl;
443 if (me == 0) std::cout <<
"Matvec original..." << std::endl;
444 origMatrix->apply(*origVector, *origProd);
447 std::cout <<
"Norm of Original matvec prod: " << origNorm << std::endl;
449 if (me == 0) std::cout <<
"Matvec redistributed..." << std::endl;
450 redistribMatrix->apply(*redistribVector, *redistribProd);
453 std::cout <<
"Norm of Redistributed matvec prod: " << redistribNorm << std::endl;
455 if (redistribNorm > origNorm+
epsilon || redistribNorm < origNorm-
epsilon) {
460 delete redistribVector;
461 delete redistribMatrix;
465 std::cout <<
"Mat-Vec product changed; FAIL" << std::endl;
469 std::cout <<
"PASS" << std::endl;
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
int main(int narg, char *arg[])
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Zoltan2::XpetraCrsGraphAdapter< SparseGraph > SparseGraphAdapter
common code used by tests
Defines the XpetraMultiVectorAdapter.
Defines XpetraCrsGraphAdapter class.
Defines the XpetraCrsMatrixAdapter class.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
Tpetra::CrsMatrix< z2TestScalar, z2TestLO, z2TestGO > SparseMatrix
InputTraits< User >::part_t part_t
An adapter for Xpetra::MultiVector.
Zoltan2::XpetraMultiVectorAdapter< IntVector > IntVectorAdapter
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
Tpetra::CrsGraph< z2TestLO, z2TestGO > SparseGraph
PartitioningProblem sets up partitioning problems for the user.
Zoltan2::XpetraCrsMatrixAdapter< SparseMatrix > SparseMatrixAdapter
Tpetra::Vector< int, z2TestLO, z2TestGO > IntVector
Defines the PartitioningProblem class.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
Zoltan2::XpetraMultiVectorAdapter< Vector > MultiVectorAdapter
void solve(bool updateInputData=true)
Direct the problem to create a solution.
std::string testDataFilePath(".")