52#include <Teuchos_ParameterList.hpp>
53#include <Teuchos_RCP.hpp>
54#include <Teuchos_CommandLineProcessor.hpp>
55#include <Tpetra_CrsMatrix.hpp>
56#include <Tpetra_Vector.hpp>
57#include <MatrixMarket_Tpetra.hpp>
79typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO>
SparseMatrix;
80typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO>
Vector;
81typedef Vector::node_type
Node;
89 typename SparseMatrix::local_inds_host_view_type indices;
90 typename SparseMatrix::values_host_view_type values;
98 for (ii=0; ii<n; ii++) {
99 A->getLocalRowView (ii, indices, values);
100 for (k=0; k< static_cast<z2TestLO>(indices.size()); k++){
104 j = perm[indices[k]];
117 return (bw_left + bw_right + 1);
132 lno_t numRows = origMatrix->getNodeNumRows();
134 std::cout <<
"origMatrix->getNodeNumRows(): " << numRows << std::endl;
137 std::cout <<
"Skipping analysis - matrix is empty" << std::endl;
141 Array<lno_t> oldMatrix(numRows*numRows);
142 Array<lno_t> newMatrix(numRows*numRows);
145 lno_t showSize = numRows;
150 std::cout << std::endl <<
"perm: ";
151 for(
lno_t n = 0; n < showSize; ++n) {
152 std::cout <<
" " << perm[n] <<
" ";
154 if(showSize < numRows) {
157 std::cout << std::endl <<
"iperm: ";
158 for(
lno_t n = 0; n < showSize; ++n) {
159 std::cout <<
" " << iperm[n] <<
" ";
161 if(showSize < numRows) {
165 std::cout << std::endl;
168 for (
lno_t j = 0; j < numRows; ++j) {
169 typename SparseMatrix::local_inds_host_view_type indices;
170 typename SparseMatrix::values_host_view_type wgts;
171 origMatrix->getLocalRowView( j, indices, wgts );
172 for (
lno_t n = 0; n < static_cast<lno_t>(indices.size()); ++n) {
173 lno_t i = indices[n];
175 oldMatrix[i + j*numRows] = 1;
176 newMatrix[perm[i] + perm[j]*numRows] = 1;
182 std::cout << std::endl <<
"original graph in matrix form:" << std::endl;
183 for(
lno_t y = 0; y < showSize; ++y) {
184 for(
lno_t x = 0; x < showSize; ++x) {
185 std::cout <<
" " << oldMatrix[x + y*numRows];
187 if(showSize < numRows) {
190 std::cout << std::endl;
192 if(showSize < numRows) {
193 for(
lno_t i = 0; i < showSize; ++i) {
198 std::cout << std::endl;
201 std::cout << std::endl <<
"reordered graph in matrix form:" << std::endl;
202 for(
lno_t y = 0; y < showSize; ++y) {
203 for(
lno_t x = 0; x < showSize; ++x) {
204 std::cout <<
" " << newMatrix[x + y*numRows];
206 if(showSize < numRows) {
209 std::cout << std::endl;
211 if(showSize < numRows) {
212 for(
lno_t i = 0; i < showSize; ++i) {
217 std::cout << std::endl;
222int mainExecute(
int narg,
char** arg, RCP<
const Teuchos::Comm<int> > comm)
224 std::string inputFile =
"";
225 std::string outputFile =
"";
226 bool verbose =
false;
229 int rank = comm->getRank();
232 Teuchos::CommandLineProcessor cmdp (
false,
false);
233 cmdp.setOption(
"inputFile", &inputFile,
234 "Name of a Matrix Market file in the data directory; "
235 "if not specified, a matrix will be generated by Galeri.");
236 cmdp.setOption(
"outputFile", &outputFile,
237 "Name of file to write the permutation");
238 cmdp.setOption(
"verbose",
"quiet", &verbose,
239 "Print messages and results.");
242 std::cout <<
"Starting everything" << std::endl;
254 std::string matrixType(
"Laplace3D");
256 cmdp.setOption(
"x", &xdim,
257 "number of gridpoints in X dimension for "
258 "mesh used to generate matrix.");
259 cmdp.setOption(
"y", &ydim,
260 "number of gridpoints in Y dimension for "
261 "mesh used to generate matrix.");
262 cmdp.setOption(
"z", &zdim,
263 "number of gridpoints in Z dimension for "
264 "mesh used to generate matrix.");
265 cmdp.setOption(
"matrix", &matrixType,
266 "Matrix type: Laplace1D, Laplace2D, or Laplace3D");
271 std::string orderMethod(
"scotch");
272 cmdp.setOption(
"order_method", &orderMethod,
273 "order_method: natural, random, rcm, scotch");
275 std::string orderMethodType(
"local");
276 cmdp.setOption(
"order_method_type", &orderMethodType,
277 "local or global or both" );
280 cmdp.parse(narg, arg);
283 RCP<UserInputForTests> uinput;
289 if (inputFile !=
""){
297 RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
300 std::cout <<
"NumRows = " << origMatrix->getGlobalNumRows() << std::endl
301 <<
"NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
302 <<
"NumProcs = " << comm->getSize() << std::endl;
317 Teuchos::ParameterList params;
318 params.set(
"order_method", orderMethod);
319 params.set(
"order_method_type", orderMethodType);
331 std::cout <<
"Going to solve" << std::endl;
342 std::cout <<
"Going to get results" << std::endl;
346 for(
int checkRank = 0; checkRank < comm->getSize(); ++checkRank ) {
348 if( checkRank == comm->getRank() ) {
349 std::cout <<
"Debug output matrix for rank: " << checkRank << std::endl;
379 if (outputFile !=
"") {
380 std::ofstream permFile;
385 std::stringstream
fname;
386 fname << outputFile <<
"." << comm->getSize() <<
"." << rank;
387 permFile.open(
fname.str().c_str());
388 for (
size_t i=0; i<checkLength; i++){
389 permFile <<
" " << checkPerm[i] << std::endl;
396 std::cout <<
"Checking permutation" << std::endl;
400 if (testReturn)
return testReturn;
404 std::cout <<
"Checking inverse permutation" << std::endl;
406 for (
size_t i=0; i< checkLength; i++){
407 testReturn = (checkInvPerm[checkPerm[i]] !=
z2TestLO(i));
408 if (testReturn)
return testReturn;
413 std::cout <<
"Checking num blocks" << std::endl;
415 testReturn = !((NumBlocks > 0) && (NumBlocks<
z2TestLO(checkLength)));
416 if (testReturn)
return testReturn;
421 std::cout <<
"Checking range" << std::endl;
423 testReturn = RangeTab[0];
424 if (testReturn)
return testReturn;
426 for (
z2TestLO i = 0; i < NumBlocks; i++){
427 testReturn = !(RangeTab[i] < RangeTab[i+1]);
428 if (testReturn)
return testReturn;
434 std::cout <<
"Checking Separator Tree" << std::endl;
438 testReturn = (TreeTab[0] != -1);
442 std::cout <<
"TreeTab[0] = " << TreeTab[0] <<
" != -1" << std::endl;
446 for (
size_t i=1; i< checkLength; i++){
447 testReturn = !(TreeTab[i] < NumBlocks);
449 std::cout <<
"TreeTab[" << i <<
"] = " << TreeTab[i] <<
" >= NumBlocks "
450 <<
" = " << NumBlocks << std::endl;
456 std::cout <<
"Going to compute the bandwidth" << std::endl;
464 std::cout <<
"Original Bandwidth: " << originalBandwidth << std::endl;
465 std::cout <<
"Permuted Bandwidth: " << permutedBandwidth << std::endl;
468 if(permutedBandwidth >= originalBandwidth) {
470 std::cout <<
"Test failed: bandwidth was not improved!" << std::endl;
472 std::cout <<
"Test in progress - returning OK for now..." << std::endl;
479 std::cout <<
"The bandwidth was improved. Good!" << std::endl;
483 catch (std::exception &e) {
484 std::cout <<
"Exception caught in ordering" << std::endl;
485 std::cout << e.what() << std::endl;
494 Tpetra::ScopeGuard tscope(&narg, &arg);
495 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
502 reduceAll<int,int>(*comm, Teuchos::EReductionType::REDUCE_MAX, 1,
503 &result, &resultReduced);
507 if (comm->getRank() == 0) {
508 std::cout <<
"Scotch Ordering test with " << comm->getSize()
509 <<
" processes has test return code " << resultReduced
510 <<
" and is exiting in the "
511 << ((resultReduced == 0 ) ?
"PASSED" :
"FAILED") <<
" state."
Defines the OrderingProblem class.
common code used by tests
Tpetra::Map ::local_ordinal_type zlno_t
std::string testDataFilePath(".")
Tpetra::Map ::global_ordinal_type zgno_t
Defines the XpetraCrsMatrixAdapter class.
InputTraits< User >::scalar_t scalar_t
InputTraits< User >::lno_t lno_t
OrderingProblem sets up ordering problems for the user.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
LocalOrderingSolution< lno_t > * getLocalOrderingSolution()
Get the local ordering solution to the problem.
index_t * getSeparatorRangeView() const
Get pointer to separator range.
bool haveSeparators() const
Do we have the separators?
index_t getNumSeparatorBlocks() const
Get number of separator column blocks.
int validatePerm()
returns 0 if permutation is valid, negative if invalid.
index_t * getPermutationView(bool inverse=false) const
Get pointer to permutation. If inverse = true, return inverse permutation. By default,...
size_t getPermutationSize() const
Get (local) size of permutation.
index_t * getSeparatorTreeView() const
Get pointer to separator tree.
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
map_t::local_ordinal_type lno_t
void tempDebugTest(RCP< SparseMatrix > origMatrix, Zoltan2::LocalOrderingSolution< z2TestLO > *soln)
int mainExecute(int narg, char **arg, RCP< const Teuchos::Comm< int > > comm)
Tpetra::CrsMatrix< z2TestScalar, z2TestLO, z2TestGO > SparseMatrix
size_t computeBandwidth(RCP< SparseMatrix > A, z2TestLO *perm)
Zoltan2::XpetraCrsMatrixAdapter< SparseMatrix > SparseMatrixAdapter
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
int main(int narg, char **arg)