Zoltan2
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 
154  try{
155  uinput =
156  rcp(new UserInputForTests(
157  testDataFilePath,std::string("simple"), comm, true));
158  }
159  catch(std::exception &e){
160  aok = false;
161  std::cout << e.what() << std::endl;
162  }
163  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
164 
165  RCP<tmatrix_t> tM; // original matrix (for checking)
166  RCP<tmatrix_t> newM; // migrated matrix
167 
168  tM = uinput->getUITpetraCrsMatrix();
169  size_t nrows = tM->getNodeNumRows();
170 
171  // To test migration in the input adapter we need a Solution object.
172 
173  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
174 
175  int nWeights = 1;
176 
179  typedef adapter_t::part_t part_t;
180 
181  part_t *p = new part_t [nrows];
182  memset(p, 0, sizeof(part_t) * nrows);
183  ArrayRCP<part_t> solnParts(p, 0, nrows, true);
184 
185  soln_t solution(env, comm, nWeights);
186  solution.setParts(solnParts);
187 
189  // User object is Tpetra::CrsMatrix
190  if (!gfail){
191  if (rank==0)
192  std::cout << "Input adapter for Tpetra::CrsMatrix" << std::endl;
193 
194  RCP<const tmatrix_t> ctM = rcp_const_cast<const tmatrix_t>(tM);
195  RCP<Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t> > tMInput;
196 
197  try {
198  tMInput =
200  }
201  catch (std::exception &e){
202  aok = false;
203  std::cout << e.what() << std::endl;
204  }
205  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter ", 1);
206 
207  fail = verifyInputAdapter<tmatrix_t>(*tMInput, *tM);
208 
209  gfail = globalFail(comm, fail);
210 
211  if (!gfail){
212  tmatrix_t *mMigrate = NULL;
213  try{
214  tMInput->applyPartitioningSolution(*tM, mMigrate, solution);
215  newM = rcp(mMigrate);
216  }
217  catch (std::exception &e){
218  fail = 11;
219  std::cout << "Error caught: " << e.what() << std::endl;
220  }
221 
222  gfail = globalFail(comm, fail);
223 
224  if (!gfail){
225  RCP<const tmatrix_t> cnewM = rcp_const_cast<const tmatrix_t>(newM);
226  RCP<Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t> > newInput;
227  try{
228  newInput = rcp(new Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t>(cnewM));
229  }
230  catch (std::exception &e){
231  aok = false;
232  std::cout << e.what() << std::endl;
233  }
234  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 2 ", 1);
235 
236  if (rank==0){
237  std::cout <<
238  "Input adapter for Tpetra::CrsMatrix migrated to proc 0" <<
239  std::endl;
240  }
241  fail = verifyInputAdapter<tmatrix_t>(*newInput, *newM);
242  if (fail) fail += 100;
243  gfail = globalFail(comm, fail);
244  }
245  }
246  if (gfail){
247  printFailureCode(comm, fail);
248  }
249  }
250 
252  // User object is Xpetra::CrsMatrix
253  if (!gfail){
254  if (rank==0)
255  std::cout << "Input adapter for Xpetra::CrsMatrix" << std::endl;
256 
257  RCP<xmatrix_t> xM = uinput->getUIXpetraCrsMatrix();
258  RCP<const xmatrix_t> cxM = rcp_const_cast<const xmatrix_t>(xM);
259  RCP<Zoltan2::XpetraCrsMatrixAdapter<xmatrix_t> > xMInput;
260 
261  try {
262  xMInput =
264  }
265  catch (std::exception &e){
266  aok = false;
267  std::cout << e.what() << std::endl;
268  }
269  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 3 ", 1);
270 
271  fail = verifyInputAdapter<xmatrix_t>(*xMInput, *tM);
272 
273  gfail = globalFail(comm, fail);
274 
275  if (!gfail){
276  xmatrix_t *mMigrate =NULL;
277  try{
278  xMInput->applyPartitioningSolution(*xM, mMigrate, solution);
279  }
280  catch (std::exception &e){
281  std::cout << "Error caught: " << e.what() << std::endl;
282  fail = 11;
283  }
284 
285  gfail = globalFail(comm, fail);
286 
287  if (!gfail){
288  RCP<const xmatrix_t> cnewM(mMigrate);
289  RCP<Zoltan2::XpetraCrsMatrixAdapter<xmatrix_t> > newInput;
290  try{
291  newInput =
293  }
294  catch (std::exception &e){
295  aok = false;
296  std::cout << e.what() << std::endl;
297  }
298  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 4 ", 1);
299 
300  if (rank==0){
301  std::cout <<
302  "Input adapter for Xpetra::CrsMatrix migrated to proc 0" <<
303  std::endl;
304  }
305  fail = verifyInputAdapter<xmatrix_t>(*newInput, *newM);
306  if (fail) fail += 100;
307  gfail = globalFail(comm, fail);
308  }
309  }
310  if (gfail){
311  printFailureCode(comm, fail);
312  }
313  }
314 
315 #ifdef HAVE_EPETRA_DATA_TYPES
316  // User object is Epetra_CrsMatrix
318  typedef Epetra_CrsMatrix ematrix_t;
319  if (!gfail){
320  if (rank==0)
321  std::cout << "Input adapter for Epetra_CrsMatrix" << std::endl;
322 
323  RCP<ematrix_t> eM = uinput->getUIEpetraCrsMatrix();
324  RCP<const ematrix_t> ceM = rcp_const_cast<const ematrix_t>(eM);
325  RCP<Zoltan2::XpetraCrsMatrixAdapter<ematrix_t> > eMInput;
326 
327  bool goodAdapter = true;
328  try {
329  eMInput =
331  }
332  catch (std::exception &e){
333  if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
334  aok = false;
335  goodAdapter = false;
336  std::cout << e.what() << std::endl;
337  }
338  else {
339  // We expect an error from Xpetra when znode_t != Xpetra::EpetraNode
340  // Ignore it, but skip tests using this matrix adapter.
341  std::cout << "Node type is not supported by Xpetra's Epetra interface;"
342  << " Skipping this test." << std::endl;
343  std::cout << "FYI: Here's the exception message: " << std::endl
344  << e.what() << std::endl;
345  goodAdapter = false;
346  }
347  }
348  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 5 ", 1);
349 
350  if (goodAdapter) {
351  fail = verifyInputAdapter<ematrix_t>(*eMInput, *tM);
352 
353  gfail = globalFail(comm, fail);
354 
355  if (!gfail){
356  ematrix_t *mMigrate =NULL;
357  try{
358  eMInput->applyPartitioningSolution(*eM, mMigrate, solution);
359  }
360  catch (std::exception &e){
361  std::cout << "Error caught: " << e.what() << std::endl;
362  fail = 11;
363  }
364 
365  gfail = globalFail(comm, fail);
366 
367  if (!gfail){
368  RCP<const ematrix_t> cnewM(mMigrate, true);
369  RCP<Zoltan2::XpetraCrsMatrixAdapter<ematrix_t> > newInput;
370  try{
371  newInput =
373  }
374  catch (std::exception &e){
375  aok = false;
376  std::cout << e.what() << std::endl;
377  }
378  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsMatrixAdapter 6 ", 1);
379 
380  if (rank==0){
381  std::cout <<
382  "Input adapter for Epetra_CrsMatrix migrated to proc 0" <<
383  std::endl;
384  }
385  fail = verifyInputAdapter<ematrix_t>(*newInput, *newM);
386  if (fail) fail += 100;
387  gfail = globalFail(comm, fail);
388  }
389  }
390  if (gfail){
391  printFailureCode(comm, fail);
392  }
393  }
394  }
395 #endif
396 
398  // DONE
399 
400  if (rank==0)
401  std::cout << "PASS" << std::endl;
402 }
Zoltan2_InputTraits.hpp
Traits for application input objects.
globalFail
int globalFail(const RCP< const Comm< int > > &comm, int fail)
Definition: ErrorHandlingForTests.hpp:109
printFailureCode
void printFailureCode(const RCP< const Comm< int > > &comm, int fail)
Definition: ErrorHandlingForTests.hpp:116
Zoltan2::PartitioningSolution
A PartitioningSolution is a solution to a partitioning problem.
Definition: Zoltan2_PartitioningSolution.hpp:55
Zoltan2::XpetraCrsMatrixAdapter::getLocalNumRows
size_t getLocalNumRows() const
Returns the number of rows on this process.
Definition: Zoltan2_XpetraCrsMatrixAdapter.hpp:162
Zoltan2_XpetraCrsMatrixAdapter.hpp
Defines the XpetraCrsMatrixAdapter class.
Zoltan2::InputTraits::offset_t
default_offset_t offset_t
The data type to represent offsets.
Definition: Zoltan2_InputTraits.hpp:195
Zoltan2::XpetraCrsMatrixAdapter
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
Definition: Zoltan2_XpetraCrsMatrixAdapter.hpp:86
testDataFilePath
std::string testDataFilePath(".")
xmatrix_t
Xpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > xmatrix_t
Definition: XpetraCrsMatrixInput.cpp:70
UserInputForTests
Definition: UserInputForTests.hpp:126
part_t
SparseMatrixAdapter_t::part_t part_t
Definition: partitioningTree.cpp:74
verifyInputAdapter
int verifyInputAdapter(Zoltan2::XpetraCrsMatrixAdapter< User > &ia, tmatrix_t &M)
Definition: XpetraCrsMatrixInput.cpp:97
Zoltan2::Environment
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
Definition: Zoltan2_Environment.hpp:83
TEST_FAIL_AND_EXIT
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Definition: ErrorHandlingForTests.hpp:70
zgno_t
int zgno_t
Definition: Zoltan2_TestHelpers.hpp:143
Zoltan2::XpetraCrsMatrixAdapter::getLocalNumColumns
size_t getLocalNumColumns() const
Returns the number of columns on this process.
Definition: Zoltan2_XpetraCrsMatrixAdapter.hpp:166
printMatrix
void printMatrix(RCP< const Comm< int > > &comm, zlno_t nrows, const zgno_t *rowIds, const offset_t *offsets, const zgno_t *colIds)
Definition: XpetraCrsMatrixInput.cpp:73
main
int main(int narg, char *arg[])
Definition: XpetraCrsMatrixInput.cpp:140
zlno_t
int zlno_t
Definition: Zoltan2_TestHelpers.hpp:142
Zoltan2::XpetraCrsMatrixAdapter::getRowIDsView
void getRowIDsView(const gno_t *&rowIds) const
Definition: Zoltan2_XpetraCrsMatrixAdapter.hpp:176
fail
static const std::string fail
Definition: findUniqueGids.cpp:80
Zoltan2_TestHelpers.hpp
common code used by tests
Zoltan2::XpetraCrsMatrixAdapter::getCRSView
void getCRSView(const offset_t *&offsets, const gno_t *&colIds) const
Definition: Zoltan2_XpetraCrsMatrixAdapter.hpp:182
tmatrix_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tmatrix_t
Definition: XpetraCrsMatrixInput.cpp:69