Zoltan2
TpetraRowMatrixInput.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::TpetraRowMatrixAdapter
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::rcp_dynamic_cast;
68 using Teuchos::Comm;
69 
70 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> ztcrsmatrix_t;
71 typedef Tpetra::RowMatrix<zscalar_t, zlno_t, zgno_t, znode_t> ztrowmatrix_t;
72 
74 
75 template<typename offset_t>
76 void printMatrix(RCP<const Comm<int> > &comm, zlno_t nrows,
77  const zgno_t *rowIds, const offset_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 (offset_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 
100 
101 template <typename User>
104 {
105  typedef typename Zoltan2::InputTraits<User>::offset_t offset_t;
106 
107  RCP<const Comm<int> > comm = M.getComm();
108  int fail = 0, gfail=0;
109 
110  if (!fail && ia.getLocalNumRows() != M.getNodeNumRows())
111  fail = 4;
112 
113  if (M.getNodeNumRows()){
114  if (!fail && ia.getLocalNumColumns() != M.getNodeNumCols())
115  fail = 6;
116  }
117 
118  gfail = globalFail(*comm, fail);
119 
120  const zgno_t *rowIds=NULL;
121  ArrayRCP<const zgno_t> colIds;
122  ArrayRCP<const offset_t> offsets;
123  size_t nrows=0;
124 
125  if (!gfail){
126 
127  nrows = ia.getLocalNumRows();
128  ia.getRowIDsView(rowIds);
129  ia.getCRSView(offsets, colIds);
130 
131  if (nrows != M.getNodeNumRows())
132  fail = 8;
133 
134  gfail = globalFail(*comm, fail);
135 
136  if (gfail == 0){
137  printMatrix<offset_t>(comm, nrows, rowIds, offsets.getRawPtr(), colIds.getRawPtr());
138  }
139  else{
140  if (!fail) fail = 10;
141  }
142  }
143  return fail;
144 }
145 
146 
147 int main(int narg, char *arg[])
148 {
149  Tpetra::ScopeGuard tscope(&narg, &arg);
150  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
151 
152  int rank = comm->getRank();
153  int fail = 0, gfail=0;
154  bool aok = true;
155 
156  // Create object that can give us Tpetra matrices for testing.
157 
158  RCP<UserInputForTests> uinput;
159  Teuchos::ParameterList params;
160  params.set("input file", "simple");
161  params.set("file type", "Chaco");
162 
163  try{
164  uinput = rcp(new UserInputForTests(params, comm));
165  }
166  catch(std::exception &e){
167  aok = false;
168  std::cout << e.what() << std::endl;
169  }
170  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
171 
172  // Input matrix and row matrix cast from it.
173  RCP<ztcrsmatrix_t> tM = uinput->getUITpetraCrsMatrix();
174  RCP<ztrowmatrix_t> trM = rcp_dynamic_cast<ztrowmatrix_t>(tM);
175 
176  RCP<ztrowmatrix_t> newM; // migrated matrix
177 
178  size_t nrows = trM->getNodeNumRows();
179 
180  // To test migration in the input adapter we need a Solution object.
181 
182  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
183 
184  int nWeights = 1;
185 
188  typedef adapter_t::part_t part_t;
189 
190  part_t *p = new part_t [nrows];
191  memset(p, 0, sizeof(part_t) * nrows);
192  ArrayRCP<part_t> solnParts(p, 0, nrows, true);
193 
194  soln_t solution(env, comm, nWeights);
195  solution.setParts(solnParts);
196 
198  // User object is Tpetra::RowMatrix
199  if (!gfail){
200  if (rank==0)
201  std::cout << "Input adapter for Tpetra::RowMatrix" << std::endl;
202 
203  RCP<const ztrowmatrix_t> ctrM = rcp_const_cast<const ztrowmatrix_t>(
204  rcp_dynamic_cast<ztrowmatrix_t>(tM));
205  RCP<adapter_t> trMInput;
206 
207  try {
208  trMInput = rcp(new adapter_t(ctrM));
209  }
210  catch (std::exception &e){
211  aok = false;
212  std::cout << e.what() << std::endl;
213  }
214  TEST_FAIL_AND_EXIT(*comm, aok, "TpetraRowMatrixAdapter ", 1);
215 
216  fail = verifyInputAdapter<ztrowmatrix_t>(*trMInput, *trM);
217 
218  gfail = globalFail(*comm, fail);
219 
220  if (!gfail){
221  ztrowmatrix_t *mMigrate = NULL;
222  try{
223  trMInput->applyPartitioningSolution(*trM, mMigrate, solution);
224  newM = rcp(mMigrate);
225  }
226  catch (std::exception &e){
227  fail = 11;
228  std::cout << "Error caught: " << e.what() << std::endl;
229  }
230 
231  gfail = globalFail(*comm, fail);
232 
233  if (!gfail){
234  RCP<const ztrowmatrix_t> cnewM =
235  rcp_const_cast<const ztrowmatrix_t>(newM);
236  RCP<adapter_t> newInput;
237  try{
238  newInput = rcp(new adapter_t(cnewM));
239  }
240  catch (std::exception &e){
241  aok = false;
242  std::cout << e.what() << std::endl;
243  }
244  TEST_FAIL_AND_EXIT(*comm, aok, "TpetraRowMatrixAdapter 2 ", 1);
245 
246  if (rank==0){
247  std::cout <<
248  "Input adapter for Tpetra::RowMatrix migrated to proc 0" <<
249  std::endl;
250  }
251  fail = verifyInputAdapter<ztrowmatrix_t>(*newInput, *newM);
252  if (fail) fail += 100;
253  gfail = globalFail(*comm, fail);
254  }
255  }
256  if (gfail){
257  printFailureCode(*comm, fail);
258  }
259  }
260 
262  // DONE
263 
264  if (rank==0)
265  std::cout << "PASS" << std::endl;
266 }
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)
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > ztcrsmatrix_t
int main(int narg, char *arg[])
int verifyInputAdapter(Zoltan2::TpetraRowMatrixAdapter< User > &ia, ztrowmatrix_t &M)
Tpetra::RowMatrix< zscalar_t, zlno_t, zgno_t, znode_t > ztrowmatrix_t
void printMatrix(RCP< const Comm< int > > &comm, zlno_t nrows, const zgno_t *rowIds, const offset_t *offsets, const zgno_t *colIds)
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 TpetraRowMatrixAdapter class.
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 Tpetra::RowMatrix data.
size_t getLocalNumColumns() const
Returns the number of columns on this process.
size_t getLocalNumRows() const
Returns the number of rows on this process.
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds) const
void getRowIDsView(const gno_t *&rowIds) const
static const std::string fail
SparseMatrixAdapter_t::part_t part_t
default_offset_t offset_t
The data type to represent offsets.