Zoltan2
Loading...
Searching...
No Matches
XpetraCrsMatrixInput.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
10//
11// Basic testing of Zoltan2::XpetraCrsMatrixAdapter
12
18#include <string>
19
23
24#include <Teuchos_DefaultComm.hpp>
25#include <Teuchos_RCP.hpp>
26#include <Teuchos_Comm.hpp>
27#include <Teuchos_CommHelpers.hpp>
28
29using Teuchos::RCP;
30using Teuchos::rcp;
31using Teuchos::rcp_const_cast;
32using Teuchos::Comm;
33
34typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> tmatrix_t;
35typedef Xpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> xmatrix_t;
36
37enum class Type {
38 CRS,
39 CCS
40};
41
42template<typename offset_t>
43void printMatrix(RCP<const Comm<int> > &comm, zlno_t nPrimaryIds,
44 const zgno_t *primaryIds, const offset_t *offsets, const zgno_t *secondaryIds, Type type)
45{
46 const auto rank = comm->getRank();
47 const auto nprocs = comm->getSize();
48 const std::string primaryIdName = type == Type::CRS ? "row" : "col";
49
50 comm->barrier();
51
52 for (int p=0; p < nprocs; p++){
53 if (p == rank){
54 std::cout << rank << ":" << std::endl;
55 for (zlno_t i=0; i < nPrimaryIds; i++){
56 std::cout << " " << primaryIdName << " " << primaryIds[i] << ": ";
57 for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
58 std::cout << secondaryIds[j] << " ";
59 }
60 std::cout << std::endl;
61 }
62 std::cout.flush();
63 }
64 comm->barrier();
65 }
66 comm->barrier();
67}
68
69template <typename User>
72{
73 typedef typename Zoltan2::InputTraits<User>::offset_t offset_t;
74
75 RCP<const Comm<int> > comm = M.getComm();
76 int fail = 0, gfail=0;
77
78 if (!fail && ia.getLocalNumRows() != M.getLocalNumRows())
79 fail = 4;
80
81 if (M.getLocalNumRows()){
82 if (!fail && ia.getLocalNumColumns() != M.getLocalNumCols())
83 fail = 6;
84 }
85
86 gfail = globalFail(*comm, fail);
87
88 const zgno_t *rowIds=NULL;
89 ArrayRCP<const zgno_t> colIds;
90 ArrayRCP<const offset_t> offsets;
91 size_t nrows=0;
92
93 if (!gfail){
94
95 nrows = ia.getLocalNumRows();
96 ia.getRowIDsView(rowIds);
97 ia.getCRSView(offsets, colIds);
98
99 if (nrows != M.getLocalNumRows())
100 fail = 8;
101
102 gfail = globalFail(*comm, fail);
103
104 if (gfail == 0){
105 printMatrix<offset_t>(comm, nrows, rowIds, offsets.getRawPtr(),
106 colIds.getRawPtr(), Type::CRS);
107 }
108 else{
109 if (!fail) fail = 10;
110 }
111
112 const auto nCols = ia.getLocalNumColumns();
113 const zgno_t *colIds2=nullptr;
114 ia.getColumnIDsView(colIds2);
115
116 ArrayRCP<const zgno_t> ccsRowIds;
117 ArrayRCP<const offset_t> ccsOffsets;
118 ia.getCCSView(ccsOffsets, ccsRowIds);
119
120 printMatrix<offset_t>(comm, nCols, colIds2, ccsOffsets.getRawPtr(),
121 ccsRowIds.getRawPtr(), Type::CCS);
122 }
123 return fail;
124}
125
126int main(int narg, char *arg[])
127{
128 Tpetra::ScopeGuard tscope(&narg, &arg);
129 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
130
131 int rank = comm->getRank();
132 int fail = 0, gfail=0;
133 bool aok = true;
134
135 // Create object that can give us test Tpetra, Xpetra
136 // and Epetra matrices for testing.
137
138 RCP<UserInputForTests> uinput;
139 Teuchos::ParameterList params;
140 params.set("input file", "simple");
141 params.set("file type", "Chaco");
142
143 try{
144 uinput = rcp(new UserInputForTests(params, comm));
145 }
146 catch(std::exception &e){
147 aok = false;
148 std::cout << e.what() << std::endl;
149 }
150 TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
151
152 RCP<tmatrix_t> tM; // original matrix (for checking)
153 RCP<tmatrix_t> newM; // migrated matrix
154
155 tM = uinput->getUITpetraCrsMatrix();
156 size_t nrows = tM->getLocalNumRows();
157
158 // To test migration in the input adapter we need a Solution object.
159
160 RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
161
162 int nWeights = 1;
163
166 typedef adapter_t::part_t part_t;
167
168 part_t *p = new part_t [nrows];
169 memset(p, 0, sizeof(part_t) * nrows);
170 ArrayRCP<part_t> solnParts(p, 0, nrows, true);
171
172 soln_t solution(env, comm, nWeights);
173 solution.setParts(solnParts);
174
176 // User object is Tpetra::CrsMatrix
177 if (!gfail){
178 if (rank==0)
179 std::cout << "Input adapter for Tpetra::CrsMatrix" << std::endl;
180
181 RCP<const tmatrix_t> ctM = rcp_const_cast<const tmatrix_t>(tM);
182 RCP<Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t> > tMInput;
183
184 try {
185 tMInput =
187 }
188 catch (std::exception &e){
189 aok = false;
190 std::cout << e.what() << std::endl;
191 }
192 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter ", 1);
193
194 fail = verifyInputAdapter<tmatrix_t>(*tMInput, *tM);
195
196 gfail = globalFail(*comm, fail);
197
198 if (!gfail){
199 tmatrix_t *mMigrate = NULL;
200 try{
201 tMInput->applyPartitioningSolution(*tM, mMigrate, solution);
202 newM = rcp(mMigrate);
203 }
204 catch (std::exception &e){
205 fail = 11;
206 std::cout << "Error caught: " << e.what() << std::endl;
207 }
208
209 gfail = globalFail(*comm, fail);
210
211 if (!gfail){
212 RCP<const tmatrix_t> cnewM = rcp_const_cast<const tmatrix_t>(newM);
213 RCP<Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t> > newInput;
214 try{
215 newInput = rcp(new Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t>(cnewM));
216 }
217 catch (std::exception &e){
218 aok = false;
219 std::cout << e.what() << std::endl;
220 }
221 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 2 ", 1);
222
223 if (rank==0){
224 std::cout <<
225 "Input adapter for Tpetra::CrsMatrix migrated to proc 0" <<
226 std::endl;
227 }
228 fail = verifyInputAdapter<tmatrix_t>(*newInput, *newM);
229 if (fail)
230 fail += 100;
231 gfail = globalFail(*comm, fail);
232 }
233 }
234 if (gfail){
235 printFailureCode(*comm, fail);
236 }
237 }
238
240 // User object is Xpetra::CrsMatrix
241 if (!gfail){
242 if (rank==0)
243 std::cout << "Input adapter for Xpetra::CrsMatrix" << std::endl;
244
245 RCP<xmatrix_t> xM = uinput->getUIXpetraCrsMatrix();
246 RCP<const xmatrix_t> cxM = rcp_const_cast<const xmatrix_t>(xM);
247 RCP<Zoltan2::XpetraCrsMatrixAdapter<xmatrix_t> > xMInput;
248
249 try {
250 xMInput =
252 }
253 catch (std::exception &e){
254 aok = false;
255 std::cout << e.what() << std::endl;
256 }
257 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 3 ", 1);
258
259 fail = verifyInputAdapter<xmatrix_t>(*xMInput, *tM);
260
261 gfail = globalFail(*comm, fail);
262
263 if (!gfail){
264 xmatrix_t *mMigrate =NULL;
265 try{
266 xMInput->applyPartitioningSolution(*xM, mMigrate, solution);
267 }
268 catch (std::exception &e){
269 std::cout << "Error caught: " << e.what() << std::endl;
270 fail = 11;
271 }
272
273 gfail = globalFail(*comm, fail);
274
275 if (!gfail){
276 RCP<const xmatrix_t> cnewM(mMigrate);
277 RCP<Zoltan2::XpetraCrsMatrixAdapter<xmatrix_t> > newInput;
278 try{
279 newInput =
281 }
282 catch (std::exception &e){
283 aok = false;
284 std::cout << e.what() << std::endl;
285 }
286 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 4 ", 1);
287
288 if (rank==0){
289 std::cout <<
290 "Input adapter for Xpetra::CrsMatrix migrated to proc 0" <<
291 std::endl;
292 }
293 fail = verifyInputAdapter<xmatrix_t>(*newInput, *newM);
294 if (fail) fail += 100;
295 gfail = globalFail(*comm, fail);
296 }
297 }
298 if (gfail){
299 printFailureCode(*comm, fail);
300 }
301 }
302
303#ifdef HAVE_EPETRA_DATA_TYPES
305 // User object is Epetra_CrsMatrix
306 typedef Epetra_CrsMatrix ematrix_t;
307 if (!gfail){
308 if (rank==0)
309 std::cout << "Input adapter for Epetra_CrsMatrix" << std::endl;
310
311 RCP<ematrix_t> eM = uinput->getUIEpetraCrsMatrix();
312 RCP<const ematrix_t> ceM = rcp_const_cast<const ematrix_t>(eM);
313 RCP<Zoltan2::XpetraCrsMatrixAdapter<ematrix_t> > eMInput;
314
315 bool goodAdapter = true;
316 try {
317 eMInput =
319 }
320 catch (std::exception &e){
321 if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
322 aok = false;
323 goodAdapter = false;
324 std::cout << e.what() << std::endl;
325 }
326 else {
327 // We expect an error from Xpetra when znode_t != Xpetra::EpetraNode
328 // Ignore it, but skip tests using this matrix adapter.
329 std::cout << "Node type is not supported by Xpetra's Epetra interface;"
330 << " Skipping this test." << std::endl;
331 std::cout << "FYI: Here's the exception message: " << std::endl
332 << e.what() << std::endl;
333 goodAdapter = false;
334 }
335 }
336 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 5 ", 1);
337
338 if (goodAdapter) {
339 fail = verifyInputAdapter<ematrix_t>(*eMInput, *tM);
340
341 gfail = globalFail(*comm, fail);
342
343 if (!gfail){
344 ematrix_t *mMigrate =NULL;
345 try{
346 eMInput->applyPartitioningSolution(*eM, mMigrate, solution);
347 }
348 catch (std::exception &e){
349 std::cout << "Error caught: " << e.what() << std::endl;
350 fail = 11;
351 }
352
353 gfail = globalFail(*comm, fail);
354
355 if (!gfail){
356 RCP<const ematrix_t> cnewM(mMigrate, true);
357 RCP<Zoltan2::XpetraCrsMatrixAdapter<ematrix_t> > newInput;
358 try{
359 newInput =
361 }
362 catch (std::exception &e){
363 aok = false;
364 std::cout << e.what() << std::endl;
365 }
366 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 6 ", 1);
367
368 if (rank==0){
369 std::cout <<
370 "Input adapter for Epetra_CrsMatrix migrated to proc 0" <<
371 std::endl;
372 }
373 fail = verifyInputAdapter<ematrix_t>(*newInput, *newM);
374 if (fail) fail += 100;
375 gfail = globalFail(*comm, fail);
376 }
377 }
378 if (gfail){
379 printFailureCode(*comm, fail);
380 }
381 }
382 }
383#endif
384
386 // DONE
387
388 if (rank==0)
389 std::cout << "PASS" << std::endl;
390}
int globalFail(const Comm< int > &comm, int fail)
void printFailureCode(const Comm< int > &comm, int fail)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Xpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > xmatrix_t
int verifyInputAdapter(Zoltan2::XpetraCrsMatrixAdapter< User > &ia, tmatrix_t &M)
void printMatrix(RCP< const Comm< int > > &comm, zlno_t nPrimaryIds, const zgno_t *primaryIds, const offset_t *offsets, const zgno_t *secondaryIds, Type type)
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tmatrix_t
Traits for application input objects.
common code used by tests
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
Defines the XpetraCrsMatrixAdapter class.
int main()
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
A PartitioningSolution is a solution to a partitioning problem.
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds) const
void getRowIDsView(const gno_t *&rowIds) const
size_t getLocalNumColumns() const
Returns the number of columns on this process.
void getCCSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &rowIds) const override
void getColumnIDsView(const gno_t *&colIds) const
size_t getLocalNumRows() const
Returns the number of rows on this process.
static const std::string fail
SparseMatrixAdapter_t::part_t part_t
default_offset_t offset_t
The data type to represent offsets.