Zoltan2
orderingScotch.cpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// Zoltan2: A package of combinatorial algorithms for scientific computing
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact Karen Devine (kddevin@sandia.gov)
39// Erik Boman (egboman@sandia.gov)
40// Siva Rajamanickam (srajama@sandia.gov)
41//
42// ***********************************************************************
43//
44// @HEADER
48#include <iostream>
49#include <fstream>
50#include <limits>
51#include <vector>
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>
58
59using Teuchos::RCP;
60
62// Program to demonstrate use of Zoltan2 to order a TPetra matrix
63// (read from a MatrixMarket file or generated by Galeri::Xpetra).
64// Usage:
65// a.out [--inputFile=filename] [--outputFile=outfile] [--verbose]
66// [--x=#] [--y=#] [--z=#] [--matrix={Laplace1D,Laplace2D,Laplace3D}
67// Karen Devine, 2011
69
71// Eventually want to use Teuchos unit tests to vary z2TestLO and
72// GO. For now, we set them at compile time based on whether Tpetra
73// is built with explicit instantiation on. (in Zoltan2_TestHelpers.hpp)
74
78
79typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO> SparseMatrix;
80typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO> Vector;
81typedef Vector::node_type Node;
83
84
85// Using perm
86size_t computeBandwidth(RCP<SparseMatrix> A, z2TestLO *perm)
87{
88 z2TestLO ii, i, j, k;
89 typename SparseMatrix::local_inds_host_view_type indices;
90 typename SparseMatrix::values_host_view_type values;
91
92 z2TestLO bw_left = 0;
93 z2TestLO bw_right = 0;
94
95 z2TestLO n = A->getNodeNumRows();
96
97 // Loop over rows of matrix
98 for (ii=0; ii<n; ii++) {
99 A->getLocalRowView (ii, indices, values);
100 for (k=0; k< static_cast<z2TestLO>(indices.size()); k++){
101 if (indices[k] < n){ // locally owned
102 if (perm){
103 i = perm[ii];
104 j = perm[indices[k]];
105 } else {
106 i = ii;
107 j = indices[k];
108 }
109 if (j-i > bw_right)
110 bw_right = j-i;
111 if (i-j > bw_left)
112 bw_left = i-j;
113 }
114 }
115 }
116 // Total bandwidth is the sum of left and right + 1
117 return (bw_left + bw_right + 1);
118}
119
120#define MDM
121#ifdef MDM
122// This is all temp code to be deleted when refactoring is compelte.
124 RCP<SparseMatrix> origMatrix, Zoltan2::LocalOrderingSolution<z2TestLO> *soln)
125{
126 typedef typename SparseMatrixAdapter::scalar_t scalar_t;
127 typedef typename SparseMatrixAdapter::lno_t lno_t;
128
129 lno_t * perm = soln->getPermutationView();
130 lno_t * iperm = soln->getPermutationView(true);
131
132 lno_t numRows = origMatrix->getNodeNumRows();
133
134 std::cout << "origMatrix->getNodeNumRows(): " << numRows << std::endl;
135
136 if (numRows == 0) {
137 std::cout << "Skipping analysis - matrix is empty" << std::endl;
138 return;
139 }
140
141 Array<lno_t> oldMatrix(numRows*numRows);
142 Array<lno_t> newMatrix(numRows*numRows);
143
144 // print the solution permutation
145 lno_t showSize = numRows;
146 if(showSize > 30) {
147 showSize = 30;
148 }
149
150 std::cout << std::endl << "perm: ";
151 for(lno_t n = 0; n < showSize; ++n) {
152 std::cout << " " << perm[n] << " ";
153 }
154 if(showSize < numRows) {
155 std::cout << " ..."; // partial print...
156 }
157 std::cout << std::endl << "iperm: ";
158 for(lno_t n = 0; n < showSize; ++n) {
159 std::cout << " " << iperm[n] << " ";
160 }
161 if(showSize < numRows) {
162 std::cout << " ..."; // partial print...
163 }
164
165 std::cout << std::endl;
166
167 // write 1's in our matrix
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];
174 if (i < numRows) { // locally owned
175 oldMatrix[i + j*numRows] = 1;
176 newMatrix[perm[i] + perm[j]*numRows] = 1;
177 }
178 }
179 }
180
181 // print oldMatrix
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];
186 }
187 if(showSize < numRows) {
188 std::cout << " ..."; // partial print...
189 }
190 std::cout << std::endl;
191 }
192 if(showSize < numRows) {
193 for(lno_t i = 0; i < showSize; ++i) {
194 std::cout << " ."; // partial print...
195 }
196 std::cout << " ..."; // partial print...
197 }
198 std::cout << std::endl;
199
200 // print newMatrix
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];
205 }
206 if(showSize < numRows) {
207 std::cout << " ..."; // partial print...
208 }
209 std::cout << std::endl;
210 }
211 if(showSize < numRows) {
212 for(lno_t i = 0; i < showSize; ++i) {
213 std::cout << " ."; // partial print...
214 }
215 std::cout << " ..."; // partial print...
216 }
217 std::cout << std::endl;
218}
219#endif
220
222int mainExecute(int narg, char** arg, RCP<const Teuchos::Comm<int> > comm)
223{
224 std::string inputFile = ""; // Matrix Market file to read
225 std::string outputFile = ""; // Output file to write
226 bool verbose = false; // Verbosity of output
227 int testReturn = 0;
228
229 int rank = comm->getRank();
230
231 // Read run-time options.
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.");
240
241 if (rank == 0 ) {
242 std::cout << "Starting everything" << std::endl;
243 }
244
246 // Even with cmdp option "true", I get errors for having these
247 // arguments on the command line. (On redsky build)
248 // KDDKDD Should just be warnings, right? Code should still work with these
249 // KDDKDD params in the create-a-matrix file. Better to have them where
250 // KDDKDD they are used.
251 int xdim=10;
252 int ydim=10;
253 int zdim=10;
254 std::string matrixType("Laplace3D");
255
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");
267
269 // Ordering options to test.
271 std::string orderMethod("scotch"); // NYI
272 cmdp.setOption("order_method", &orderMethod,
273 "order_method: natural, random, rcm, scotch");
274
275 std::string orderMethodType("local");
276 cmdp.setOption("order_method_type", &orderMethodType,
277 "local or global or both" );
278
280 cmdp.parse(narg, arg);
281
282
283 RCP<UserInputForTests> uinput;
284
285 // MDM - temp testing code
286 // testDataFilePath = "/Users/micheldemessieres/Documents/trilinos-build/trilinosrepo/parZoltan2/packages/zoltan2/test/driver/driverinputs/ordering";
287 // inputFile = "orderingTest";
288
289 if (inputFile != ""){ // Input file specified; read a matrix
290 uinput = rcp(new UserInputForTests(testDataFilePath, inputFile, comm, true));
291 }
292 else { // Let Galeri generate a matrix
293 uinput = rcp(
294 new UserInputForTests(xdim, ydim, zdim, matrixType, comm, true, true));
295 }
296
297 RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
298
299 if (rank == 0) {
300 std::cout << "NumRows = " << origMatrix->getGlobalNumRows() << std::endl
301 << "NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
302 << "NumProcs = " << comm->getSize() << std::endl;
303 }
304
306 // Currently Not Used
307 /*
308 RCP<Vector> origVector, origProd;
309 origProd = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
310 origMatrix->getRangeMap());
311 origVector = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
312 origMatrix->getDomainMap());
313 origVector->randomize();
314 */
315
317 Teuchos::ParameterList params;
318 params.set("order_method", orderMethod);
319 params.set("order_method_type", orderMethodType);
320
322 SparseMatrixAdapter adapter(origMatrix);
323
325
326 try {
327 Zoltan2::OrderingProblem<SparseMatrixAdapter> problem(&adapter, &params,
328 comm);
329
330 if( rank == 0 ) {
331 std::cout << "Going to solve" << std::endl;
332 }
333 problem.solve();
334
336 size_t checkLength;
337 z2TestLO *checkPerm, *checkInvPerm;
339 problem.getLocalOrderingSolution();
340
341 if( rank == 0 ) {
342 std::cout << "Going to get results" << std::endl;
343 }
344
345 #ifdef MDM // Temp debugging code all of which can be removed for final
346 for( int checkRank = 0; checkRank < comm->getSize(); ++checkRank ) {
347 comm->barrier();
348 if( checkRank == comm->getRank() ) {
349 std::cout << "Debug output matrix for rank: " << checkRank << std::endl;
350 tempDebugTest(origMatrix, soln);
351 }
352 comm->barrier();
353 }
354 #endif
355
356 // Permutation
357 checkLength = soln->getPermutationSize();
358 checkPerm = soln->getPermutationView();
359 checkInvPerm = soln->getPermutationView(true); // get permutation inverse
360
361 // Separators.
362 // The following methods needs to be supported:
363 // haveSeparators: true if Scotch Nested Dissection was called.
364 // getCBlkPtr: *CBlkPtr from Scotch_graphOrder
365 // getRangeTab: RangeTab from Scotch_graphOrder
366 // getTreeTab: TreeTab from Scotch_graphOrder
367 z2TestLO NumBlocks = 0;
368 z2TestLO *RangeTab;
369 z2TestLO *TreeTab;
370 if (soln->haveSeparators()) {
371 NumBlocks = soln->getNumSeparatorBlocks(); // BDD
372 RangeTab = soln->getSeparatorRangeView(); // BDD
373 TreeTab = soln->getSeparatorTreeView(); // BDD
374 }
375 else {
376 // TODO FAIL with error
377 }
378
379 if (outputFile != "") {
380 std::ofstream permFile;
381
382 // Write permutation (0-based) to file
383 // each process writes local perm to a separate file
384 //std::string fname = outputFile + "." + me;
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;
390 }
391 permFile.close();
392 }
393
394 // Validate that checkPerm is a permutation
395 if (rank == 0 ) {
396 std::cout << "Checking permutation" << std::endl;
397 }
398
399 testReturn = soln->validatePerm();
400 if (testReturn) return testReturn;
401
402 // Validate the inverse permutation.
403 if (rank == 0 ) {
404 std::cout << "Checking inverse permutation" << std::endl;
405 }
406 for (size_t i=0; i< checkLength; i++){
407 testReturn = (checkInvPerm[checkPerm[i]] != z2TestLO(i));
408 if (testReturn) return testReturn;
409 }
410
411 // Validate NumBlocks
412 if (rank == 0 ) {
413 std::cout << "Checking num blocks" << std::endl;
414 }
415 testReturn = !((NumBlocks > 0) && (NumBlocks<z2TestLO(checkLength)));
416 if (testReturn) return testReturn;
417
418 // Validate RangeTab.
419 // Should be monitonically increasing, RT[0] = 0; RT[NumBlocks+1]=nVtx;
420 if (rank == 0 ) {
421 std::cout << "Checking range" << std::endl;
422 }
423 testReturn = RangeTab[0];
424 if (testReturn) return testReturn;
425
426 for (z2TestLO i = 0; i < NumBlocks; i++){
427 testReturn = !(RangeTab[i] < RangeTab[i+1]);
428 if (testReturn) return testReturn;
429 }
430
431 // How do we validate TreeTab?
432 // TreeTab root has -1, other values < NumBlocks
433 if (rank == 0 ) {
434 std::cout << "Checking Separator Tree" << std::endl;
435 }
436
437 if (checkLength) {
438 testReturn = (TreeTab[0] != -1);
439 }
440
441 if (testReturn) {
442 std::cout << "TreeTab[0] = " << TreeTab[0] << " != -1" << std::endl;
443 return testReturn;
444 }
445
446 for (size_t i=1; i< checkLength; i++){
447 testReturn = !(TreeTab[i] < NumBlocks);
448 if (testReturn) {
449 std::cout << "TreeTab[" << i << "] = " << TreeTab[i] << " >= NumBlocks "
450 << " = " << NumBlocks << std::endl;
451 return testReturn;
452 }
453 }
454
455 if (rank == 0 ) {
456 std::cout << "Going to compute the bandwidth" << std::endl;
457 }
458
459 // Compute original and permuted bandwidth
460 z2TestLO originalBandwidth = computeBandwidth(origMatrix, nullptr);
461 z2TestLO permutedBandwidth = computeBandwidth(origMatrix, checkPerm);
462
463 if (rank == 0 ) {
464 std::cout << "Original Bandwidth: " << originalBandwidth << std::endl;
465 std::cout << "Permuted Bandwidth: " << permutedBandwidth << std::endl;
466 }
467
468 if(permutedBandwidth >= originalBandwidth) {
469 if (rank == 0 ) {
470 std::cout << "Test failed: bandwidth was not improved!" << std::endl;
471
472 std::cout << "Test in progress - returning OK for now..." << std::endl;
473 }
474
475 // return 1; // TO DO - we need test to have proper ordering
476 }
477 else {
478 if (rank == 0) {
479 std::cout << "The bandwidth was improved. Good!" << std::endl;
480 }
481 }
482 }
483 catch (std::exception &e) {
484 std::cout << "Exception caught in ordering" << std::endl;
485 std::cout << e.what() << std::endl;
486 return 1;
487 }
488
489 return 0;
490}
491
492int main(int narg, char** arg)
493{
494 Tpetra::ScopeGuard tscope(&narg, &arg);
495 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
496
497 int result = mainExecute(narg, arg, comm);
498
499 // get reduced return code form all procsses
500 comm->barrier();
501 int resultReduced;
502 reduceAll<int,int>(*comm, Teuchos::EReductionType::REDUCE_MAX, 1,
503 &result, &resultReduced);
504
505 // provide a final message which guarantees that the test will fail
506 // if any of the processes failed
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."
512 << std::endl;
513 }
514}
515
Defines the OrderingProblem class.
common code used by tests
float zscalar_t
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.
zlno_t z2TestLO
Definition: coloring1.cpp:75
map_t::local_ordinal_type lno_t
Definition: mapRemotes.cpp:17
fname
Begin.
Definition: validXML.py:19
Vector::node_type Node
zlno_t z2TestLO
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
zgno_t z2TestGO
zscalar_t z2TestScalar
int main(int narg, char **arg)