Zoltan2
Loading...
Searching...
No Matches
orderingScotch.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
13#include <iostream>
14#include <fstream>
15#include <limits>
16#include <vector>
17#include <Teuchos_ParameterList.hpp>
18#include <Teuchos_RCP.hpp>
19#include <Teuchos_CommandLineProcessor.hpp>
20#include <Tpetra_CrsMatrix.hpp>
21#include <Tpetra_Vector.hpp>
22#include <MatrixMarket_Tpetra.hpp>
23
24using Teuchos::RCP;
25
27// Program to demonstrate use of Zoltan2 to order a TPetra matrix
28// (read from a MatrixMarket file or generated by Galeri::Xpetra).
29// Usage:
30// a.out [--inputFile=filename] [--outputFile=outfile] [--verbose]
31// [--x=#] [--y=#] [--z=#] [--matrix={Laplace1D,Laplace2D,Laplace3D}
32// Karen Devine, 2011
34
36// Eventually want to use Teuchos unit tests to vary z2TestLO and
37// GO. For now, we set them at compile time based on whether Tpetra
38// is built with explicit instantiation on. (in Zoltan2_TestHelpers.hpp)
39
43
44typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO> SparseMatrix;
45typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO> Vector;
46typedef Vector::node_type Node;
48
49
50// Using perm
51size_t computeBandwidth(RCP<SparseMatrix> A, z2TestLO *perm)
52{
53 z2TestLO ii, i, j, k;
54 typename SparseMatrix::local_inds_host_view_type indices;
55 typename SparseMatrix::values_host_view_type values;
56
57 z2TestLO bw_left = 0;
58 z2TestLO bw_right = 0;
59
61
62 // Loop over rows of matrix
63 for (ii=0; ii<n; ii++) {
64 A->getLocalRowView (ii, indices, values);
65 for (k=0; k< static_cast<z2TestLO>(indices.size()); k++){
66 if (indices[k] < n){ // locally owned
67 if (perm){
68 i = perm[ii];
69 j = perm[indices[k]];
70 } else {
71 i = ii;
72 j = indices[k];
73 }
74 if (j-i > bw_right)
75 bw_right = j-i;
76 if (i-j > bw_left)
77 bw_left = i-j;
78 }
79 }
80 }
81 // Total bandwidth is the sum of left and right + 1
82 return (bw_left + bw_right + 1);
83}
84
85#define MDM
86#ifdef MDM
87// This is all temp code to be deleted when refactoring is compelte.
89 RCP<SparseMatrix> origMatrix, Zoltan2::LocalOrderingSolution<z2TestLO> *soln)
90{
91 typedef typename SparseMatrixAdapter::lno_t lno_t;
92
93 lno_t * perm = soln->getPermutationView();
94 lno_t * iperm = soln->getPermutationView(true);
95
96 lno_t numRows = origMatrix->getLocalNumRows();
97
98 std::cout << "origMatrix->getLocalNumRows(): " << numRows << std::endl;
99
100 if (numRows == 0) {
101 std::cout << "Skipping analysis - matrix is empty" << std::endl;
102 return;
103 }
104
105 Array<lno_t> oldMatrix(numRows*numRows);
106 Array<lno_t> newMatrix(numRows*numRows);
107
108 // print the solution permutation
109 lno_t showSize = numRows;
110 if(showSize > 30) {
111 showSize = 30;
112 }
113
114 std::cout << std::endl << "perm: ";
115 for(lno_t n = 0; n < showSize; ++n) {
116 std::cout << " " << perm[n] << " ";
117 }
118 if(showSize < numRows) {
119 std::cout << " ..."; // partial print...
120 }
121 std::cout << std::endl << "iperm: ";
122 for(lno_t n = 0; n < showSize; ++n) {
123 std::cout << " " << iperm[n] << " ";
124 }
125 if(showSize < numRows) {
126 std::cout << " ..."; // partial print...
127 }
128
129 std::cout << std::endl;
130
131 // write 1's in our matrix
132 for (lno_t j = 0; j < numRows; ++j) {
133 typename SparseMatrix::local_inds_host_view_type indices;
134 typename SparseMatrix::values_host_view_type wgts;
135 origMatrix->getLocalRowView( j, indices, wgts );
136 for (lno_t n = 0; n < static_cast<lno_t>(indices.size()); ++n) {
137 lno_t i = indices[n];
138 if (i < numRows) { // locally owned
139 oldMatrix[i + j*numRows] = 1;
140 newMatrix[perm[i] + perm[j]*numRows] = 1;
141 }
142 }
143 }
144
145 // print oldMatrix
146 std::cout << std::endl << "original graph in matrix form:" << std::endl;
147 for(lno_t y = 0; y < showSize; ++y) {
148 for(lno_t x = 0; x < showSize; ++x) {
149 std::cout << " " << oldMatrix[x + y*numRows];
150 }
151 if(showSize < numRows) {
152 std::cout << " ..."; // partial print...
153 }
154 std::cout << std::endl;
155 }
156 if(showSize < numRows) {
157 for(lno_t i = 0; i < showSize; ++i) {
158 std::cout << " ."; // partial print...
159 }
160 std::cout << " ..."; // partial print...
161 }
162 std::cout << std::endl;
163
164 // print newMatrix
165 std::cout << std::endl << "reordered graph in matrix form:" << std::endl;
166 for(lno_t y = 0; y < showSize; ++y) {
167 for(lno_t x = 0; x < showSize; ++x) {
168 std::cout << " " << newMatrix[x + y*numRows];
169 }
170 if(showSize < numRows) {
171 std::cout << " ..."; // partial print...
172 }
173 std::cout << std::endl;
174 }
175 if(showSize < numRows) {
176 for(lno_t i = 0; i < showSize; ++i) {
177 std::cout << " ."; // partial print...
178 }
179 std::cout << " ..."; // partial print...
180 }
181 std::cout << std::endl;
182}
183#endif
184
186int mainExecute(int narg, char** arg, RCP<const Teuchos::Comm<int> > comm)
187{
188 std::string inputFile = ""; // Matrix Market file to read
189 std::string outputFile = ""; // Output file to write
190 bool verbose = false; // Verbosity of output
191 int testReturn = 0;
192
193 int rank = comm->getRank();
194
195 // Read run-time options.
196 Teuchos::CommandLineProcessor cmdp (false, false);
197 cmdp.setOption("inputFile", &inputFile,
198 "Name of a Matrix Market file in the data directory; "
199 "if not specified, a matrix will be generated by Galeri.");
200 cmdp.setOption("outputFile", &outputFile,
201 "Name of file to write the permutation");
202 cmdp.setOption("verbose", "quiet", &verbose,
203 "Print messages and results.");
204
205 if (rank == 0 ) {
206 std::cout << "Starting everything" << std::endl;
207 }
208
210 // Even with cmdp option "true", I get errors for having these
211 // arguments on the command line. (On redsky build)
212 // KDDKDD Should just be warnings, right? Code should still work with these
213 // KDDKDD params in the create-a-matrix file. Better to have them where
214 // KDDKDD they are used.
215 int xdim=10;
216 int ydim=10;
217 int zdim=10;
218 std::string matrixType("Laplace3D");
219
220 cmdp.setOption("x", &xdim,
221 "number of gridpoints in X dimension for "
222 "mesh used to generate matrix.");
223 cmdp.setOption("y", &ydim,
224 "number of gridpoints in Y dimension for "
225 "mesh used to generate matrix.");
226 cmdp.setOption("z", &zdim,
227 "number of gridpoints in Z dimension for "
228 "mesh used to generate matrix.");
229 cmdp.setOption("matrix", &matrixType,
230 "Matrix type: Laplace1D, Laplace2D, or Laplace3D");
231
233 // Ordering options to test.
235 std::string orderMethod("scotch"); // NYI
236 cmdp.setOption("order_method", &orderMethod,
237 "order_method: natural, random, rcm, scotch");
238
239 std::string orderMethodType("local");
240 cmdp.setOption("order_method_type", &orderMethodType,
241 "local or global or both" );
242
244 cmdp.parse(narg, arg);
245
246
247 RCP<UserInputForTests> uinput;
248
249 // MDM - temp testing code
250 // testDataFilePath = "/Users/micheldemessieres/Documents/trilinos-build/trilinosrepo/parZoltan2/packages/zoltan2/test/driver/driverinputs/ordering";
251 // inputFile = "orderingTest";
252
253 if (inputFile != ""){ // Input file specified; read a matrix
254 uinput = rcp(new UserInputForTests(testDataFilePath, inputFile, comm, true));
255 }
256 else { // Let Galeri generate a matrix
257 uinput = rcp(
258 new UserInputForTests(xdim, ydim, zdim, matrixType, comm, true, true));
259 }
260
261 RCP<SparseMatrix> origMatrix = uinput->getUITpetraCrsMatrix();
262
263 if (rank == 0) {
264 std::cout << "NumRows = " << origMatrix->getGlobalNumRows() << std::endl
265 << "NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
266 << "NumProcs = " << comm->getSize() << std::endl;
267 }
268
270 // Currently Not Used
271 /*
272 RCP<Vector> origVector, origProd;
273 origProd = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
274 origMatrix->getRangeMap());
275 origVector = Tpetra::createVector<z2TestScalar,z2TestLO,z2TestGO>(
276 origMatrix->getDomainMap());
277 origVector->randomize();
278 */
279
281 Teuchos::ParameterList params;
282 params.set("order_method", orderMethod);
283 params.set("order_method_type", orderMethodType);
284
286 SparseMatrixAdapter adapter(origMatrix);
287
289
290 try {
291 Zoltan2::OrderingProblem<SparseMatrixAdapter> problem(&adapter, &params,
292 comm);
293
294 if( rank == 0 ) {
295 std::cout << "Going to solve" << std::endl;
296 }
297 problem.solve();
298
300 size_t checkLength;
301 z2TestLO *checkPerm, *checkInvPerm;
303 problem.getLocalOrderingSolution();
304
305 if( rank == 0 ) {
306 std::cout << "Going to get results" << std::endl;
307 }
308
309 #ifdef MDM // Temp debugging code all of which can be removed for final
310 for( int checkRank = 0; checkRank < comm->getSize(); ++checkRank ) {
311 comm->barrier();
312 if( checkRank == comm->getRank() ) {
313 std::cout << "Debug output matrix for rank: " << checkRank << std::endl;
314 tempDebugTest(origMatrix, soln);
315 }
316 comm->barrier();
317 }
318 #endif
319
320 // Permutation
321 checkLength = soln->getPermutationSize();
322 checkPerm = soln->getPermutationView();
323 checkInvPerm = soln->getPermutationView(true); // get permutation inverse
324
325 // Separators.
326 // The following methods needs to be supported:
327 // haveSeparators: true if Scotch Nested Dissection was called.
328 // getCBlkPtr: *CBlkPtr from Scotch_graphOrder
329 // getRangeTab: RangeTab from Scotch_graphOrder
330 // getTreeTab: TreeTab from Scotch_graphOrder
331 z2TestLO NumBlocks = 0;
332 z2TestLO *RangeTab;
333 z2TestLO *TreeTab;
334 if (soln->haveSeparators()) {
335 NumBlocks = soln->getNumSeparatorBlocks(); // BDD
336 RangeTab = soln->getSeparatorRangeView(); // BDD
337 TreeTab = soln->getSeparatorTreeView(); // BDD
338 }
339 else {
340 // TODO FAIL with error
341 }
342
343 if (outputFile != "") {
344 std::ofstream permFile;
345
346 // Write permutation (0-based) to file
347 // each process writes local perm to a separate file
348 //std::string fname = outputFile + "." + me;
349 std::stringstream fname;
350 fname << outputFile << "." << comm->getSize() << "." << rank;
351 permFile.open(fname.str().c_str());
352 for (size_t i=0; i<checkLength; i++){
353 permFile << " " << checkPerm[i] << std::endl;
354 }
355 permFile.close();
356 }
357
358 // Validate that checkPerm is a permutation
359 if (rank == 0 ) {
360 std::cout << "Checking permutation" << std::endl;
361 }
362
363 testReturn = soln->validatePerm();
364 if (testReturn) return testReturn;
365
366 // Validate the inverse permutation.
367 if (rank == 0 ) {
368 std::cout << "Checking inverse permutation" << std::endl;
369 }
370 for (size_t i=0; i< checkLength; i++){
371 testReturn = (checkInvPerm[checkPerm[i]] != z2TestLO(i));
372 if (testReturn) return testReturn;
373 }
374
375 // Validate NumBlocks
376 if (rank == 0 ) {
377 std::cout << "Checking num blocks" << std::endl;
378 }
379 testReturn = !((NumBlocks > 0) && (NumBlocks<z2TestLO(checkLength)));
380 if (testReturn) return testReturn;
381
382 // Validate RangeTab.
383 // Should be monitonically increasing, RT[0] = 0; RT[NumBlocks+1]=nVtx;
384 if (rank == 0 ) {
385 std::cout << "Checking range" << std::endl;
386 }
387 testReturn = RangeTab[0];
388 if (testReturn) return testReturn;
389
390 for (z2TestLO i = 0; i < NumBlocks; i++){
391 testReturn = !(RangeTab[i] < RangeTab[i+1]);
392 if (testReturn) return testReturn;
393 }
394
395 // How do we validate TreeTab?
396 // TreeTab root has -1, other values < NumBlocks
397 if (rank == 0 ) {
398 std::cout << "Checking Separator Tree" << std::endl;
399 }
400
401 if (checkLength) {
402 testReturn = (TreeTab[0] != -1);
403 }
404
405 if (testReturn) {
406 std::cout << "TreeTab[0] = " << TreeTab[0] << " != -1" << std::endl;
407 return testReturn;
408 }
409
410 for (size_t i=1; i< checkLength; i++){
411 testReturn = !(TreeTab[i] < NumBlocks);
412 if (testReturn) {
413 std::cout << "TreeTab[" << i << "] = " << TreeTab[i] << " >= NumBlocks "
414 << " = " << NumBlocks << std::endl;
415 return testReturn;
416 }
417 }
418
419 if (rank == 0 ) {
420 std::cout << "Going to compute the bandwidth" << std::endl;
421 }
422
423 // Compute original and permuted bandwidth
424 z2TestLO originalBandwidth = computeBandwidth(origMatrix, nullptr);
425 z2TestLO permutedBandwidth = computeBandwidth(origMatrix, checkPerm);
426
427 if (rank == 0 ) {
428 std::cout << "Original Bandwidth: " << originalBandwidth << std::endl;
429 std::cout << "Permuted Bandwidth: " << permutedBandwidth << std::endl;
430 }
431
432 if(permutedBandwidth >= originalBandwidth) {
433 if (rank == 0 ) {
434 std::cout << "Test failed: bandwidth was not improved!" << std::endl;
435
436 std::cout << "Test in progress - returning OK for now..." << std::endl;
437 }
438
439 // return 1; // TO DO - we need test to have proper ordering
440 }
441 else {
442 if (rank == 0) {
443 std::cout << "The bandwidth was improved. Good!" << std::endl;
444 }
445 }
446 }
447 catch (std::exception &e) {
448 std::cout << "Exception caught in ordering" << std::endl;
449 std::cout << e.what() << std::endl;
450 return 1;
451 }
452
453 return 0;
454}
455
456int main(int narg, char** arg)
457{
458 Tpetra::ScopeGuard tscope(&narg, &arg);
459 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
460
461 int result = mainExecute(narg, arg, comm);
462
463 // get reduced return code form all procsses
464 comm->barrier();
465 int resultReduced;
466 reduceAll<int,int>(*comm, Teuchos::EReductionType::REDUCE_MAX, 1,
467 &result, &resultReduced);
468
469 // provide a final message which guarantees that the test will fail
470 // if any of the processes failed
471 if (comm->getRank() == 0) {
472 std::cout << "Scotch Ordering test with " << comm->getSize()
473 << " processes has test return code " << resultReduced
474 << " and is exiting in the "
475 << ((resultReduced == 0 ) ? "PASSED" : "FAILED") << " state."
476 << std::endl;
477 }
478}
479
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.
int main()
typename 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.
size_t getLocalNumRows() const
Returns the number of rows on this process.
zlno_t z2TestLO
Definition coloring1.cpp:40
map_t::local_ordinal_type lno_t
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