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_GlobalMPISession.hpp>
60 #include <Teuchos_DefaultComm.hpp>
61 #include <Teuchos_RCP.hpp>
62 #include <Teuchos_Comm.hpp>
63 #include <Teuchos_CommHelpers.hpp>
64 
65 using namespace std;
66 using Teuchos::RCP;
67 using Teuchos::rcp;
68 using Teuchos::rcp_const_cast;
69 using Teuchos::Comm;
70 using Teuchos::DefaultComm;
71 
72 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> tmatrix_t;
73 typedef Xpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> xmatrix_t;
74 typedef Epetra_CrsMatrix ematrix_t;
75 
76 void printMatrix(RCP<const Comm<int> > &comm, zlno_t nrows,
77  const zgno_t *rowIds, const zlno_t *offsets, const zgno_t *colIds)
78 {
79  int rank = comm->getRank();
80  int nprocs = comm->getSize();
81  comm->barrier();
82  for (int p=0; p < nprocs; p++){
83  if (p == rank){
84  std::cout << rank << ":" << std::endl;
85  for (zlno_t i=0; i < nrows; i++){
86  std::cout << " row " << rowIds[i] << ": ";
87  for (zlno_t j=offsets[i]; j < offsets[i+1]; j++){
88  std::cout << colIds[j] << " ";
89  }
90  std::cout << std::endl;
91  }
92  std::cout.flush();
93  }
94  comm->barrier();
95  }
96  comm->barrier();
97 }
98 
99 template <typename User>
102 {
103  RCP<const Comm<int> > comm = M.getComm();
104  int fail = 0, gfail=0;
105 
106  if (!fail && ia.getLocalNumRows() != M.getNodeNumRows())
107  fail = 4;
108 
109  if (M.getNodeNumRows()){
110  if (!fail && ia.getLocalNumColumns() != M.getNodeNumCols())
111  fail = 6;
112  }
113 
114  gfail = globalFail(comm, fail);
115 
116  const zgno_t *rowIds=NULL, *colIds=NULL;
117  const zlno_t *offsets=NULL;
118  size_t nrows=0;
119 
120  if (!gfail){
121 
122  nrows = ia.getLocalNumRows();
123  ia.getRowIDsView(rowIds);
124  ia.getCRSView(offsets, colIds);
125 
126  if (nrows != M.getNodeNumRows())
127  fail = 8;
128 
129  gfail = globalFail(comm, fail);
130 
131  if (gfail == 0){
132  printMatrix(comm, nrows, rowIds, offsets, colIds);
133  }
134  else{
135  if (!fail) fail = 10;
136  }
137  }
138  return fail;
139 }
140 
141 int main(int argc, char *argv[])
142 {
143  Teuchos::GlobalMPISession session(&argc, &argv);
144  RCP<const Comm<int> > comm = DefaultComm<int>::getComm();
145  int rank = comm->getRank();
146  int fail = 0, gfail=0;
147 
148  // Create object that can give us test Tpetra, Xpetra
149  // and Epetra matrices for testing.
150 
151  RCP<UserInputForTests> uinput;
152 
153  try{
154  uinput =
155  rcp(new UserInputForTests(
156  testDataFilePath,std::string("simple"), comm, true));
157  }
158  catch(std::exception &e){
159  TEST_FAIL_AND_EXIT(*comm, 0, string("input ")+e.what(), 1);
160  }
161 
162  RCP<tmatrix_t> tM; // original matrix (for checking)
163  RCP<tmatrix_t> newM; // migrated matrix
164 
165  tM = uinput->getUITpetraCrsMatrix();
166  size_t nrows = tM->getNodeNumRows();
167 
168  // To test migration in the input adapter we need a Solution object.
169 
170  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment);
171 
172  int nWeights = 1;
173 
176  typedef adapter_t::part_t part_t;
177 
178  part_t *p = new part_t [nrows];
179  memset(p, 0, sizeof(part_t) * nrows);
180  ArrayRCP<part_t> solnParts(p, 0, nrows, true);
181 
182  soln_t solution(env, comm, nWeights);
183  solution.setParts(solnParts);
184 
186  // User object is Tpetra::CrsMatrix
187  if (!gfail){
188  RCP<const tmatrix_t> ctM = rcp_const_cast<const tmatrix_t>(tM);
189  RCP<Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t> > tMInput;
190 
191  try {
192  tMInput =
194  }
195  catch (std::exception &e){
196  TEST_FAIL_AND_EXIT(*comm, 0,
197  string("XpetraCrsMatrixAdapter ")+e.what(), 1);
198  }
199 
200  if (rank==0)
201  std::cout << "Input adapter for Tpetra::CrsMatrix" << std::endl;
202 
203  fail = verifyInputAdapter<tmatrix_t>(*tMInput, *tM);
204 
205  gfail = globalFail(comm, fail);
206 
207  if (!gfail){
208  tmatrix_t *mMigrate = NULL;
209  try{
210  tMInput->applyPartitioningSolution(*tM, mMigrate, solution);
211  newM = rcp(mMigrate);
212  }
213  catch (std::exception &e){
214  fail = 11;
215  cout << "Error caught: " << e.what() << endl;
216  }
217 
218  gfail = globalFail(comm, fail);
219 
220  if (!gfail){
221  RCP<const tmatrix_t> cnewM = rcp_const_cast<const tmatrix_t>(newM);
222  RCP<Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t> > newInput;
223  try{
224  newInput = rcp(new Zoltan2::XpetraCrsMatrixAdapter<tmatrix_t>(cnewM));
225  }
226  catch (std::exception &e){
227  TEST_FAIL_AND_EXIT(*comm, 0,
228  string("XpetraCrsMatrixAdapter 2 ")+e.what(), 1);
229  }
230 
231  if (rank==0){
232  std::cout <<
233  "Input adapter for Tpetra::CrsMatrix migrated to proc 0" <<
234  std::endl;
235  }
236  fail = verifyInputAdapter<tmatrix_t>(*newInput, *newM);
237  if (fail) fail += 100;
238  gfail = globalFail(comm, fail);
239  }
240  }
241  if (gfail){
242  printFailureCode(comm, fail);
243  }
244  }
245 
247  // User object is Xpetra::CrsMatrix
248  if (!gfail){
249  RCP<xmatrix_t> xM = uinput->getUIXpetraCrsMatrix();
250  RCP<const xmatrix_t> cxM = rcp_const_cast<const xmatrix_t>(xM);
251  RCP<Zoltan2::XpetraCrsMatrixAdapter<xmatrix_t> > xMInput;
252 
253  try {
254  xMInput =
256  }
257  catch (std::exception &e){
258  TEST_FAIL_AND_EXIT(*comm, 0,
259  string("XpetraCrsMatrixAdapter 3 ")+e.what(), 1);
260  }
261 
262  if (rank==0){
263  std::cout << "Input adapter for Xpetra::CrsMatrix" << std::endl;
264  }
265  fail = verifyInputAdapter<xmatrix_t>(*xMInput, *tM);
266 
267  gfail = globalFail(comm, fail);
268 
269  if (!gfail){
270  xmatrix_t *mMigrate =NULL;
271  try{
272  xMInput->applyPartitioningSolution(*xM, mMigrate, solution);
273  }
274  catch (std::exception &e){
275  cout << "Error caught: " << e.what() << endl;
276  fail = 11;
277  }
278 
279  gfail = globalFail(comm, fail);
280 
281  if (!gfail){
282  RCP<const xmatrix_t> cnewM(mMigrate);
283  RCP<Zoltan2::XpetraCrsMatrixAdapter<xmatrix_t> > newInput;
284  try{
285  newInput =
287  }
288  catch (std::exception &e){
289  TEST_FAIL_AND_EXIT(*comm, 0,
290  string("XpetraCrsMatrixAdapter 4 ")+e.what(), 1);
291  }
292 
293  if (rank==0){
294  std::cout <<
295  "Input adapter for Xpetra::CrsMatrix migrated to proc 0" <<
296  std::endl;
297  }
298  fail = verifyInputAdapter<xmatrix_t>(*newInput, *newM);
299  if (fail) fail += 100;
300  gfail = globalFail(comm, fail);
301  }
302  }
303  if (gfail){
304  printFailureCode(comm, fail);
305  }
306  }
307 
308 #ifdef HAVE_EPETRA_DATA_TYPES
309  // User object is Epetra_CrsMatrix
311  if (!gfail){
312  RCP<ematrix_t> eM = uinput->getUIEpetraCrsMatrix();
313  RCP<const ematrix_t> ceM = rcp_const_cast<const ematrix_t>(eM);
314  RCP<Zoltan2::XpetraCrsMatrixAdapter<ematrix_t> > eMInput;
315 
316  try {
317  eMInput =
319  }
320  catch (std::exception &e){
321  TEST_FAIL_AND_EXIT(*comm, 0,
322  string("XpetraCrsMatrixAdapter 5 ")+e.what(), 1);
323  }
324 
325  if (rank==0){
326  std::cout << "Input adapter for Epetra_CrsMatrix" << std::endl;
327  }
328  fail = verifyInputAdapter<ematrix_t>(*eMInput, *tM);
329 
330  gfail = globalFail(comm, fail);
331 
332  if (!gfail){
333  ematrix_t *mMigrate =NULL;
334  try{
335  eMInput->applyPartitioningSolution(*eM, mMigrate, solution);
336  }
337  catch (std::exception &e){
338  cout << "Error caught: " << e.what() << endl;
339  fail = 11;
340  }
341 
342  gfail = globalFail(comm, fail);
343 
344  if (!gfail){
345  RCP<const ematrix_t> cnewM(mMigrate, true);
346  RCP<Zoltan2::XpetraCrsMatrixAdapter<ematrix_t> > newInput;
347  try{
348  newInput =
350  }
351  catch (std::exception &e){
352  TEST_FAIL_AND_EXIT(*comm, 0,
353  string("XpetraCrsMatrixAdapter 6 ")+e.what(), 1);
354  }
355 
356  if (rank==0){
357  std::cout <<
358  "Input adapter for Epetra_CrsMatrix migrated to proc 0" <<
359  std::endl;
360  }
361  fail = verifyInputAdapter<ematrix_t>(*newInput, *newM);
362  if (fail) fail += 100;
363  gfail = globalFail(comm, fail);
364  }
365  }
366  if (gfail){
367  printFailureCode(comm, fail);
368  }
369  }
370 #endif
371 
373  // DONE
374 
375  if (rank==0)
376  std::cout << "PASS" << std::endl;
377 }
int globalFail(const RCP< 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
Xpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > xmatrix_t
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
void getCRSView(const lno_t *&offsets, const gno_t *&colIds) const
int main(int argc, char *argv[])
size_t getLocalNumColumns() const
Returns the number of columns on this process.
int zlno_t
common code used by tests
Defines the XpetraCrsMatrixAdapter class.
void printMatrix(RCP< const Comm< int > > &comm, zlno_t nrows, const zgno_t *rowIds, const zlno_t *offsets, const zgno_t *colIds)
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...
int zgno_t
int verifyInputAdapter(Zoltan2::XpetraCrsMatrixAdapter< User > &ia, tmatrix_t &M)
void printFailureCode(const RCP< const Comm< int > > &comm, int fail)
Epetra_CrsMatrix ematrix_t
void getRowIDsView(const gno_t *&rowIds) const
std::string testDataFilePath(".")