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