Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
XpetraCrsGraphInput.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Zoltan2: A package of combinatorial algorithms for scientific computing
4 //
5 // Copyright 2012 NTESS and the Zoltan2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 //
11 // Basic testing of Zoltan2::XpetraCrsGraphAdapter
17 #include <string>
18 
21 #include <Zoltan2_TestHelpers.hpp>
22 
23 #include <Teuchos_DefaultComm.hpp>
24 #include <Teuchos_RCP.hpp>
25 #include <Teuchos_Comm.hpp>
26 #include <Teuchos_CommHelpers.hpp>
27 
28 using Teuchos::RCP;
29 using Teuchos::rcp;
30 using Teuchos::rcp_const_cast;
31 using Teuchos::Comm;
32 using Teuchos::Array;
33 using Teuchos::ArrayView;
34 
35 typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t> tgraph_t;
36 typedef Xpetra::CrsGraph<zlno_t, zgno_t, znode_t> xgraph_t;
37 
38 template<typename offset_t>
39 void printGraph(RCP<const Comm<int> > &comm, zlno_t nvtx,
40  const zgno_t *vtxIds, const offset_t *offsets, const zgno_t *edgeIds)
41 {
42  int rank = comm->getRank();
43  int nprocs = comm->getSize();
44  comm->barrier();
45  for (int p=0; p < nprocs; p++){
46  if (p == rank){
47  std::cout << rank << ":" << std::endl;
48  for (zlno_t i=0; i < nvtx; i++){
49  std::cout << " vertex " << vtxIds[i] << ": ";
50  for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
51  std::cout << edgeIds[j] << " ";
52  }
53  std::cout << std::endl;
54  }
55  std::cout.flush();
56  }
57  comm->barrier();
58  }
59  comm->barrier();
60 }
61 
62 template <typename User>
65  tgraph_t &graph
66 )
67 {
68  typedef typename Zoltan2::InputTraits<User>::offset_t offset_t;
69 
70  RCP<const Comm<int> > comm = graph.getComm();
71  int fail = 0, gfail=0;
72 
73  if (!fail &&
74  ia.getLocalNumVertices() != graph.getLocalNumRows())
75  fail = 4;
76 
77  if (!fail &&
78  ia.getLocalNumEdges() != graph.getLocalNumEntries())
79  fail = 6;
80 
81  gfail = globalFail(*comm, fail);
82 
83  const zgno_t *vtxIds=NULL, *edgeIds=NULL;
84  const offset_t *offsets=NULL;
85  size_t nvtx=0;
86 
87  if (!gfail){
88 
89  nvtx = ia.getLocalNumVertices();
90  ia.getVertexIDsView(vtxIds);
91  ia.getEdgesView(offsets, edgeIds);
92 
93  if (nvtx != graph.getLocalNumRows())
94  fail = 8;
95 
96  gfail = globalFail(*comm, fail);
97 
98  if (gfail == 0){
99  printGraph<offset_t>(comm, nvtx, vtxIds, offsets, edgeIds);
100  }
101  else{
102  if (!fail) fail = 10;
103  }
104  }
105  return fail;
106 }
107 
108 int main(int narg, char *arg[])
109 {
110  Tpetra::ScopeGuard tscope(&narg, &arg);
111  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
112 
113  int rank = comm->getRank();
114  int fail = 0, gfail=0;
115  bool aok = true;
116 
117  // Create an object that can give us test Tpetra, Xpetra
118  // and Epetra graphs for testing.
119 
120  RCP<UserInputForTests> uinput;
121  Teuchos::ParameterList params;
122  params.set("input file", "simple");
123  params.set("file type", "Chaco");
124 
125  try{
126  uinput = rcp(new UserInputForTests(params, comm));
127  }
128  catch(std::exception &e){
129  aok = false;
130  std::cout << e.what() << std::endl;
131  }
132  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
133 
134  RCP<tgraph_t> tG; // original graph (for checking)
135  RCP<tgraph_t> newG; // migrated graph
136 
137  tG = uinput->getUITpetraCrsGraph();
138  size_t nvtx = tG->getLocalNumRows();
139 
140  // To test migration in the input adapter we need a Solution object.
141  // Our solution just assigns all objects to part zero.
142 
143  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
144 
145  int nWeights = 1;
146 
149  typedef adapter_t::part_t part_t;
150 
151  part_t *p = new part_t [nvtx];
152  memset(p, 0, sizeof(part_t) * nvtx);
153  ArrayRCP<part_t> solnParts(p, 0, nvtx, true);
154 
155  soln_t solution(env, comm, nWeights);
156  solution.setParts(solnParts);
157 
159  // User object is Tpetra::CrsGraph
160  if (!gfail){
161  if (rank==0)
162  std::cout << "Input adapter for Tpetra::CrsGraph" << std::endl;
163 
164  RCP<const tgraph_t> ctG = rcp_const_cast<const tgraph_t>(tG);
165  RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > tGInput;
166 
167  try {
168  tGInput =
170  }
171  catch (std::exception &e){
172  aok = false;
173  std::cout << e.what() << std::endl;
174  }
175  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter ", 1);
176 
177  fail = verifyInputAdapter<tgraph_t>(*tGInput, *tG);
178 
179  gfail = globalFail(*comm, fail);
180 
181  if (!gfail){
182  tgraph_t *mMigrate = NULL;
183  try{
184  tGInput->applyPartitioningSolution( *tG, mMigrate, solution);
185  newG = rcp(mMigrate);
186  }
187  catch (std::exception &e){
188  fail = 11;
189  }
190 
191  gfail = globalFail(*comm, fail);
192 
193  if (!gfail){
194  RCP<const tgraph_t> cnewG = rcp_const_cast<const tgraph_t>(newG);
195  RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > newInput;
196  try{
197  newInput = rcp(new Zoltan2::XpetraCrsGraphAdapter<tgraph_t>(cnewG));
198  }
199  catch (std::exception &e){
200  aok = false;
201  std::cout << e.what() << std::endl;
202  }
203  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 2 ", 1);
204 
205  if (rank==0){
206  std::cout <<
207  "Input adapter for Tpetra::CrsGraph migrated to proc 0" <<
208  std::endl;
209  }
210  fail = verifyInputAdapter<tgraph_t>(*newInput, *newG);
211  if (fail) fail += 100;
212  gfail = globalFail(*comm, fail);
213  }
214  }
215  if (gfail){
216  printFailureCode(*comm, fail);
217  }
218  }
219 
221  // User object is Xpetra::CrsGraph
222  if (!gfail){
223  if (rank==0)
224  std::cout << "Input adapter for Xpetra::CrsGraph" << std::endl;
225 
226  RCP<xgraph_t> xG = uinput->getUIXpetraCrsGraph();
227  RCP<const xgraph_t> cxG = rcp_const_cast<const xgraph_t>(xG);
228  RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > xGInput;
229 
230  try {
231  xGInput =
233  }
234  catch (std::exception &e){
235  aok = false;
236  std::cout << e.what() << std::endl;
237  }
238  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 3 ", 1);
239 
240  fail = verifyInputAdapter<xgraph_t>(*xGInput, *tG);
241 
242  gfail = globalFail(*comm, fail);
243 
244  if (!gfail){
245  xgraph_t *mMigrate =NULL;
246  try{
247  xGInput->applyPartitioningSolution( *xG, mMigrate, solution);
248  }
249  catch (std::exception &e){
250  fail = 11;
251  }
252 
253  gfail = globalFail(*comm, fail);
254 
255  if (!gfail){
256  RCP<const xgraph_t> cnewG(mMigrate);
257  RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > newInput;
258  try{
259  newInput =
261  }
262  catch (std::exception &e){
263  aok = false;
264  std::cout << e.what() << std::endl;
265  }
266  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 4 ", 1);
267 
268  if (rank==0){
269  std::cout <<
270  "Input adapter for Xpetra::CrsGraph migrated to proc 0" <<
271  std::endl;
272  }
273  fail = verifyInputAdapter<xgraph_t>(*newInput, *newG);
274  if (fail) fail += 100;
275  gfail = globalFail(*comm, fail);
276  }
277  }
278  if (gfail){
279  printFailureCode(*comm, fail);
280  }
281  }
282 
283 #ifdef HAVE_EPETRA_DATA_TYPES
284  // User object is Epetra_CrsGraph
286  typedef Epetra_CrsGraph egraph_t;
287  if (!gfail){
288  if (rank==0)
289  std::cout << "Input adapter for Epetra_CrsGraph" << std::endl;
290 
291  RCP<egraph_t> eG = uinput->getUIEpetraCrsGraph();
292  RCP<const egraph_t> ceG = rcp_const_cast<const egraph_t>(eG);
293  RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > eGInput;
294 
295  bool goodAdapter = true;
296  try {
297  eGInput =
299  }
300  catch (std::exception &e){
301  if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
302  aok = false;
303  goodAdapter = false;
304  std::cout << e.what() << std::endl;
305  }
306  else {
307  // We expect an error from Xpetra when znode_t != Xpetra::EpetraNode
308  // Ignore it, but skip tests using this matrix adapter.
309  std::cout << "Node type is not supported by Xpetra's Epetra interface;"
310  << " Skipping this test." << std::endl;
311  std::cout << "FYI: Here's the exception message: " << std::endl
312  << e.what() << std::endl;
313  goodAdapter = false;
314  }
315  }
316  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 5 ", 1);
317 
318  if (goodAdapter) {
319  fail = verifyInputAdapter<egraph_t>(*eGInput, *tG);
320 
321  gfail = globalFail(*comm, fail);
322 
323  if (!gfail){
324  egraph_t *mMigrate =NULL;
325  try{
326  eGInput->applyPartitioningSolution( *eG, mMigrate, solution);
327  }
328  catch (std::exception &e){
329  fail = 11;
330  }
331 
332  gfail = globalFail(*comm, fail);
333 
334  if (!gfail){
335  RCP<const egraph_t> cnewG(mMigrate, true);
336  RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > newInput;
337  try{
338  newInput =
340  }
341  catch (std::exception &e){
342  aok = false;
343  std::cout << e.what() << std::endl;
344  }
345  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 6 ", 1);
346 
347  if (rank==0){
348  std::cout <<
349  "Input adapter for Epetra_CrsGraph migrated to proc 0" <<
350  std::endl;
351  }
352  fail = verifyInputAdapter<egraph_t>(*newInput, *newG);
353  if (fail) fail += 100;
354  gfail = globalFail(*comm, fail);
355  }
356  }
357  if (gfail){
358  printFailureCode(*comm, fail);
359  }
360  }
361  }
362 #endif
363 
365  // DONE
366 
367  if (rank==0)
368  std::cout << "PASS" << std::endl;
369 }
void printFailureCode(const Comm< int > &comm, int fail)
void verifyInputAdapter(adapter_t &ia, matrix_t &matrix)
default_offset_t offset_t
The data type to represent offsets.
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Xpetra::CrsGraph< zlno_t, zgno_t, znode_t > xgraph_t
int main(int narg, char **arg)
Definition: coloring1.cpp:164
size_t getLocalNumVertices() const override
Returns the number of vertices on this process.
Defines the PartitioningSolution class.
common code used by tests
void printGraph(RCP< const Comm< int > > &comm, zlno_t nvtx, const zgno_t *vtxIds, const offset_t *offsets, const zgno_t *edgeIds)
Defines XpetraCrsGraphAdapter class.
SparseMatrixAdapter_t::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
size_t getLocalNumEdges() const override
Returns the number of edges on this process.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
Tpetra::Map::local_ordinal_type zlno_t
static const std::string fail
int globalFail(const Comm< int > &comm, int fail)
void getVertexIDsView(const gno_t *&ids) const override
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const override
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tgraph_t
Tpetra::Map::global_ordinal_type zgno_t