Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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, *colIds=NULL;
121  const offset_t *offsets=NULL;
122  size_t nrows=0;
123 
124  if (!gfail){
125 
126  nrows = ia.getLocalNumRows();
127  ia.getRowIDsView(rowIds);
128  ia.getCRSView(offsets, colIds);
129 
130  if (nrows != M.getNodeNumRows())
131  fail = 8;
132 
133  gfail = globalFail(*comm, fail);
134 
135  if (gfail == 0){
136  printMatrix<offset_t>(comm, nrows, rowIds, offsets, colIds);
137  }
138  else{
139  if (!fail) fail = 10;
140  }
141  }
142  return fail;
143 }
144 
145 
146 int main(int narg, char *arg[])
147 {
148  Tpetra::ScopeGuard tscope(&narg, &arg);
149  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
150 
151  int rank = comm->getRank();
152  int fail = 0, gfail=0;
153  bool aok = true;
154 
155  // Create object that can give us Tpetra matrices for testing.
156 
157  RCP<UserInputForTests> uinput;
158  Teuchos::ParameterList params;
159  params.set("input file", "simple");
160  params.set("file type", "Chaco");
161 
162  try{
163  uinput = rcp(new UserInputForTests(params, comm));
164  }
165  catch(std::exception &e){
166  aok = false;
167  std::cout << e.what() << std::endl;
168  }
169  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
170 
171  // Input matrix and row matrix cast from it.
172  RCP<ztcrsmatrix_t> tM = uinput->getUITpetraCrsMatrix();
173  RCP<ztrowmatrix_t> trM = rcp_dynamic_cast<ztrowmatrix_t>(tM);
174 
175  RCP<ztrowmatrix_t> newM; // migrated matrix
176 
177  size_t nrows = trM->getNodeNumRows();
178 
179  // To test migration in the input adapter we need a Solution object.
180 
181  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
182 
183  int nWeights = 1;
184 
187  typedef adapter_t::part_t part_t;
188 
189  part_t *p = new part_t [nrows];
190  memset(p, 0, sizeof(part_t) * nrows);
191  ArrayRCP<part_t> solnParts(p, 0, nrows, true);
192 
193  soln_t solution(env, comm, nWeights);
194  solution.setParts(solnParts);
195 
197  // User object is Tpetra::RowMatrix
198  if (!gfail){
199  if (rank==0)
200  std::cout << "Input adapter for Tpetra::RowMatrix" << std::endl;
201 
202  RCP<const ztrowmatrix_t> ctrM = rcp_const_cast<const ztrowmatrix_t>(
203  rcp_dynamic_cast<ztrowmatrix_t>(tM));
204  RCP<adapter_t> trMInput;
205 
206  try {
207  trMInput = rcp(new adapter_t(ctrM));
208  }
209  catch (std::exception &e){
210  aok = false;
211  std::cout << e.what() << std::endl;
212  }
213  TEST_FAIL_AND_EXIT(*comm, aok, "TpetraRowMatrixAdapter ", 1);
214 
215  fail = verifyInputAdapter<ztrowmatrix_t>(*trMInput, *trM);
216 
217  gfail = globalFail(*comm, fail);
218 
219  if (!gfail){
220  ztrowmatrix_t *mMigrate = NULL;
221  try{
222  trMInput->applyPartitioningSolution(*trM, mMigrate, solution);
223  newM = rcp(mMigrate);
224  }
225  catch (std::exception &e){
226  fail = 11;
227  std::cout << "Error caught: " << e.what() << std::endl;
228  }
229 
230  gfail = globalFail(*comm, fail);
231 
232  if (!gfail){
233  RCP<const ztrowmatrix_t> cnewM =
234  rcp_const_cast<const ztrowmatrix_t>(newM);
235  RCP<adapter_t> newInput;
236  try{
237  newInput = rcp(new adapter_t(cnewM));
238  }
239  catch (std::exception &e){
240  aok = false;
241  std::cout << e.what() << std::endl;
242  }
243  TEST_FAIL_AND_EXIT(*comm, aok, "TpetraRowMatrixAdapter 2 ", 1);
244 
245  if (rank==0){
246  std::cout <<
247  "Input adapter for Tpetra::RowMatrix migrated to proc 0" <<
248  std::endl;
249  }
250  fail = verifyInputAdapter<ztrowmatrix_t>(*newInput, *newM);
251  if (fail) fail += 100;
252  gfail = globalFail(*comm, fail);
253  }
254  }
255  if (gfail){
256  printFailureCode(*comm, fail);
257  }
258  }
259 
261  // DONE
262 
263  if (rank==0)
264  std::cout << "PASS" << std::endl;
265 }
void printFailureCode(const Comm< int > &comm, int fail)
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > ztcrsmatrix_t
default_offset_t offset_t
The data type to represent offsets.
int main(int narg, char *arg[])
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
int zlno_t
common code used by tests
Provides access for Zoltan2 to Tpetra::RowMatrix data.
SparseMatrixAdapter_t::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
void getRowIDsView(const gno_t *&rowIds) const
Sets pointer to this process&#39; rows&#39; global IDs.
size_t getLocalNumColumns() const
Returns the number of columns on this process.
Traits for application input objects.
Defines the TpetraRowMatrixAdapter class.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
void getCRSView(const offset_t *&offsets, const gno_t *&colIds) const
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)
Tpetra::RowMatrix< zscalar_t, zlno_t, zgno_t, znode_t > ztrowmatrix_t
size_t getLocalNumRows() const
Returns the number of rows on this process.
void printMatrix(RCP< const Comm< int > > &comm, zlno_t nrows, const zgno_t *rowIds, const offset_t *offsets, const zgno_t *colIds)