Zoltan2
Loading...
Searching...
No Matches
nd.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
14#include <iostream>
15#include <limits>
16#include <Teuchos_ParameterList.hpp>
17#include <Teuchos_RCP.hpp>
18#include <Teuchos_FancyOStream.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// Eventually want to use Teuchos unit tests to vary z2TestLO and
28// GO. For now, we set them at compile time based on whether Tpetra
29// is built with explicit instantiation on. (in Zoltan2_TestHelpers.hpp)
30
34
35typedef Tpetra::CrsMatrix<z2TestScalar, z2TestLO, z2TestGO> SparseMatrix_t;
36typedef Tpetra::Vector<z2TestScalar, z2TestLO, z2TestGO> Vector;
37typedef Vector::node_type Node;
38
39typedef Tpetra::MultiVector<z2TestScalar, z2TestLO, z2TestGO,znode_t> tMVector_t;
40
41
42//typedef Zoltan2::XpetraCrsMatrixAdapter<SparseMatrix_t> SparseMatrixAdapter;
44
46
47#define epsilon 0.00000001
48
49int testNDwithRCB(RCP<SparseMatrix_t> &origMatrix,RCP<tMVector_t> &coords, int numParts, int me);
50int testNDwithPHG(RCP<SparseMatrix_t> &origMatrix,int numParts, int me);
51
54int main(int narg, char** arg)
55{
59 Tpetra::ScopeGuard tscope(&narg, &arg);
60 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
61 int me = comm->getRank();
63
64 std::string inputFile = ""; // Matrix Market or Zoltan file to read
65 std::string inputPath = testDataFilePath; // Directory with input file
66 bool distributeInput = true;
67 int success = 0;
68 int numParts = 2;
69
71 // Read run-time options.
73 Teuchos::CommandLineProcessor cmdp (false, false);
74 cmdp.setOption("inputPath", &inputPath,
75 "Path to the MatrixMarket or Zoltan file to be read; "
76 "if not specified, a default path will be used.");
77 cmdp.setOption("inputFile", &inputFile,
78 "Name of the Matrix Market or Zoltan file to read; "
79 "");
80 cmdp.setOption("distribute", "no-distribute", &distributeInput,
81 "for Zoltan input files only, "
82 "indicate whether or not to distribute "
83 "input across the communicator");
84 cmdp.setOption("numParts", &numParts,
85 "Global number of parts;");
86
87 Teuchos::CommandLineProcessor::EParseCommandLineReturn
88 parseReturn= cmdp.parse( narg, arg );
89
90 if( parseReturn == Teuchos::CommandLineProcessor::PARSE_HELP_PRINTED )
91 {
92 return 0;
93 }
95
97 // Construct matrix from file
99 RCP<UserInputForTests> uinput;
100
101 if (inputFile != "") // Input file specified; read a matrix
102 {
103 uinput = rcp(new UserInputForTests(inputPath, inputFile, comm,
104 true, distributeInput));
105 }
106 else
107 {
108 std::cout << "Input file must be specified." << std::endl;
109 }
110
111 RCP<SparseMatrix_t> origMatrix = uinput->getUITpetraCrsMatrix();
112
113
114 if (me == 0)
115 {
116 std::cout << "NumRows = " << origMatrix->getGlobalNumRows() << std::endl
117 << "NumNonzeros = " << origMatrix->getGlobalNumEntries() << std::endl
118 << "NumProcs = " << comm->getSize() << std::endl
119 << "NumParts = " << numParts << std::endl;
120 }
121
122 if (origMatrix->getGlobalNumRows() < 40)
123 {
124 Teuchos::FancyOStream out(Teuchos::rcp(&std::cout,false));
125 origMatrix->describe(out, Teuchos::VERB_EXTREME);
126 }
128
130 // Get coordinates, corresponding to graph vertices
132 RCP<tMVector_t> coords;
133 try
134 {
135 coords = uinput->getUICoordinates();
136 }
137 catch(...)
138 {
139 if (me == 0)
140 std::cout << "FAIL: get coordinates" << std::endl;
141 return 1;
142 }
144
146 // Test ND ordering with RCB to compute the separator
148 testNDwithRCB(origMatrix,coords,numParts,me);
150
152 // Test ND ordering with PHG to compute the separator
154 testNDwithPHG(origMatrix,numParts,me);
156
157
158
159 return success;
160}
162
165int testNDwithRCB(RCP<SparseMatrix_t> &origMatrix,RCP<tMVector_t> &coords, int numParts, int me)
166{
167 int success=0;
168
172 Teuchos::ParameterList params;
173
174 params.set("num_global_parts", numParts);
175 params.set("order_method", "nd");
176 params.set("edge_separator_method", "rcb");
178
182 SparseMatrixAdapter_t matAdapter(origMatrix);
183
184 MultiVectorAdapter_t *ca = NULL;
185
186 try
187 {
188 ca = new MultiVectorAdapter_t(coords);
189 }
190 catch(...)
191 {
192 if (me == 0)
193 std::cout << "FAIL: vector adapter" << std::endl;
194 return 1;
195 }
196
197 matAdapter.setCoordinateInput(ca);
199
203 Zoltan2::OrderingProblem<SparseMatrixAdapter_t> problem(&matAdapter, &params);
204
205
206 try
207 {
208 if (me == 0) std::cout << "Calling solve() " << std::endl;
209 problem.solve();
210 if (me == 0) std::cout << "Done solve() " << std::endl;
211 }
212 catch (std::runtime_error &e)
213 {
214 std::cout << "Runtime exception returned from solve(): " << e.what();
215 if (!strncmp(e.what(), "BUILD ERROR", 11)) {
216 // Catching build errors as exceptions is OK in the tests
217 std::cout << " PASS" << std::endl;
218 return 0;
219 }
220 else {
221 // All other runtime_errors are failures
222 std::cout << " FAIL" << std::endl;
223 return -1;
224 }
225 }
226 catch (std::logic_error &e)
227 {
228 std::cout << "Logic exception returned from solve(): " << e.what()
229 << " FAIL" << std::endl;
230 return -1;
231 }
232 catch (std::bad_alloc &e)
233 {
234 std::cout << "Bad_alloc exception returned from solve(): " << e.what()
235 << " FAIL" << std::endl;
236 return -1;
237 }
238 catch (std::exception &e)
239 {
240 std::cout << "Unknown exception returned from solve(). " << e.what()
241 << " FAIL" << std::endl;
242 return -1;
243 }
244
246
247 delete ca;
248
249
250 std::cout << "PASS" << std::endl;
251 return success;
252}
254
255
258int testNDwithPHG(RCP<SparseMatrix_t> &origMatrix,int numParts, int me)
259{
260 int success=0;
261
265 Teuchos::ParameterList params;
266
267 params.set("num_global_parts", numParts);
268 params.set("order_method", "nd");
269 params.set("edge_separator_method", "phg");
271
275 SparseMatrixAdapter_t matAdapter(origMatrix);
277
281 Zoltan2::OrderingProblem<SparseMatrixAdapter_t> problem(&matAdapter, &params);
282
283
284 try
285 {
286 if (me == 0) std::cout << "Calling solve() " << std::endl;
287 problem.solve();
288 if (me == 0) std::cout << "Done solve() " << std::endl;
289 }
290 catch (std::runtime_error &e)
291 {
292 std::cout << "Runtime exception returned from solve(): " << e.what();
293 if (!strncmp(e.what(), "BUILD ERROR", 11)) {
294 // Catching build errors as exceptions is OK in the tests
295 std::cout << " PASS" << std::endl;
296 return 0;
297 }
298 else {
299 // All other runtime_errors are failures
300 std::cout << " FAIL" << std::endl;
301 return -1;
302 }
303 }
304 catch (std::logic_error &e)
305 {
306 std::cout << "Logic exception returned from solve(): " << e.what()
307 << " FAIL" << std::endl;
308 return -1;
309 }
310 catch (std::bad_alloc &e)
311 {
312 std::cout << "Bad_alloc exception returned from solve(): " << e.what()
313 << " FAIL" << std::endl;
314 return -1;
315 }
316 catch (std::exception &e)
317 {
318 std::cout << "Unknown exception returned from solve(). " << e.what()
319 << " FAIL" << std::endl;
320 return -1;
321 }
323
324 std::cout << "PASS" << std::endl;
325 return success;
326}
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.
Defines the XpetraMultiVectorAdapter.
int main()
void setCoordinateInput(VectorAdapter< UserCoord > *coordData) override
Allow user to provide additional data that contains coordinate info associated with the MatrixAdapter...
OrderingProblem sets up ordering problems for the user.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tMVector_t
Vector::node_type Node
Definition nd.cpp:37
int testNDwithPHG(RCP< SparseMatrix_t > &origMatrix, int numParts, int me)
Definition nd.cpp:258
zlno_t z2TestLO
Definition nd.cpp:31
Tpetra::CrsMatrix< z2TestScalar, z2TestLO, z2TestGO > SparseMatrix_t
Definition nd.cpp:35
int testNDwithRCB(RCP< SparseMatrix_t > &origMatrix, RCP< tMVector_t > &coords, int numParts, int me)
Definition nd.cpp:165
Zoltan2::XpetraMultiVectorAdapter< tMVector_t > MultiVectorAdapter_t
Definition nd.cpp:45
Zoltan2::XpetraCrsMatrixAdapter< SparseMatrix_t, tMVector_t > SparseMatrixAdapter_t
Definition nd.cpp:43
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
Definition nd.cpp:36
zgno_t z2TestGO
Definition nd.cpp:32
zscalar_t z2TestScalar
Definition nd.cpp:33