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