50 #include <Teuchos_ParameterList.hpp> 51 #include <Teuchos_RCP.hpp> 52 #include <Teuchos_CommandLineProcessor.hpp> 53 #include <Tpetra_DefaultPlatform.hpp> 54 #include <Tpetra_CrsMatrix.hpp> 55 #include <Tpetra_Vector.hpp> 56 #include <Galeri_XpetraMaps.hpp> 57 #include <Galeri_XpetraProblemFactory.hpp> 67 int main(
int narg,
char** arg)
70 Teuchos::GlobalMPISession mpiSession(&narg, &arg, NULL);
71 RCP<const Teuchos::Comm<int> > comm =
72 Tpetra::DefaultPlatform::getDefaultPlatform().getComm();
73 int me = comm->getRank();
81 typedef Tpetra::Map<zlno_t, zgno_t> Map_t;
82 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t> Matrix_t;
83 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t> MultiVector_t;
84 typedef Tpetra::Vector<zscalar_t, zlno_t, zgno_t> Vector_t;
91 std::string method =
"scotch";
92 zgno_t nx = 50, ny = 40, nz = 30;
96 Teuchos::CommandLineProcessor cmdp (
false,
false);
97 cmdp.setOption(
"method", &method,
98 "Partitioning method to use: scotch or parmetis.");
99 cmdp.setOption(
"nx", &nx,
100 "number of gridpoints in X dimension for " 101 "mesh used to generate matrix; must be >= 1.");
102 cmdp.setOption(
"ny", &ny,
103 "number of gridpoints in Y dimension for " 104 "mesh used to generate matrix; must be >= 1.");
105 cmdp.setOption(
"nz", &nz,
106 "number of gridpoints in Z dimension for " 107 "mesh used to generate matrix; must be >= 1.");
108 cmdp.parse(narg, arg);
110 if ((nx < 1) || (ny < 1) || (nz < 1)) {
111 cout <<
"Input error: nx, ny and nz must be >= 1" << endl;
119 Teuchos::ParameterList galeriList;
120 galeriList.set(
"nx", nx);
121 galeriList.set(
"ny", ny);
122 galeriList.set(
"nz", nz);
125 RCP<Matrix_t> origMatrix;
128 RCP<const Map_t> map = rcp(
new Map_t(nGlobalElements, 0, comm));
130 typedef Galeri::Xpetra::Problem<Map_t,Matrix_t,MultiVector_t> Galeri_t;
131 RCP<Galeri_t> galeriProblem =
133 Map_t, Matrix_t, MultiVector_t>
134 (
"Laplace3D", map, galeriList);
135 origMatrix = galeriProblem->BuildMatrix();
137 catch (std::exception &e) {
138 cout <<
"Exception in Galeri matrix generation. " << e.what() << endl;
143 cout <<
"NumRows = " << origMatrix->getGlobalNumRows() << endl
144 <<
"NumNonzeros = " << origMatrix->getGlobalNumEntries() << endl
145 <<
"NumProcs = " << comm->getSize() << endl;
148 RCP<Vector_t> origVector, origProd;
149 origProd = Tpetra::createVector<zscalar_t,zlno_t,zgno_t>(
150 origMatrix->getRangeMap());
151 origVector = Tpetra::createVector<zscalar_t,zlno_t,zgno_t>(
152 origMatrix->getDomainMap());
153 origVector->randomize();
156 Teuchos::ParameterList param;
157 param.set(
"partitioning_approach",
"partition");
158 param.set(
"algorithm", method);
161 MatrixAdapter_t adapter(origMatrix);
169 catch (std::exception &e) {
170 cout <<
"Exception returned from solve(). " << e.what() << endl;
177 if (me == 0) cout <<
"Redistributing matrix..." << endl;
178 RCP<Matrix_t> redistribMatrix;
179 adapter.applyPartitioningSolution(*origMatrix, redistribMatrix,
182 if (me == 0) cout <<
"Redistributing vectors..." << endl;
183 RCP<Vector_t> redistribVector;
184 MultiVectorAdapter_t adapterVector(origVector);
185 adapterVector.applyPartitioningSolution(*origVector, redistribVector,
189 RCP<Vector_t> redistribProd;
190 redistribProd = Tpetra::createVector<zscalar_t,zlno_t,zgno_t>(
191 redistribMatrix->getRangeMap());
195 if (origMatrix->getGlobalNumRows() <= 50) {
196 cout << me <<
" ORIGINAL: ";
197 for (
size_t i = 0; i < origVector->getLocalLength(); i++)
198 cout << origVector->getMap()->getGlobalElement(i) <<
" ";
200 cout << me <<
" REDISTRIB: ";
201 for (
size_t i = 0; i < redistribVector->getLocalLength(); i++)
202 cout << redistribVector->getMap()->getGlobalElement(i) <<
" ";
210 if (me == 0) cout <<
"Matvec original..." << endl;
211 origMatrix->apply(*origVector, *origProd);
212 zscalar_t origNorm = origProd->norm2();
214 cout <<
"Norm of Original matvec prod: " << origNorm << endl;
216 if (me == 0) cout <<
"Matvec redistributed..." << endl;
217 redistribMatrix->apply(*redistribVector, *redistribProd);
218 zscalar_t redistribNorm = redistribProd->norm2();
220 cout <<
"Norm of Redistributed matvec prod: " << redistribNorm << endl;
223 const double epsilon = 0.00000001;
224 if (redistribNorm > origNorm+epsilon || redistribNorm < origNorm-epsilon)
225 cout <<
"Mat-Vec product changed; FAIL" << std::endl;
227 cout <<
"PASS" << endl;
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
Defines the XpetraMultiVectorAdapter.
Defines the XpetraCrsMatrixAdapter class.
int main(int narg, char **arg)
An adapter for Xpetra::MultiVector.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
void solve(bool updateInputData=true)
Direct the problem to create a solution.