Zoltan2
XpetraCrsGraphInput.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::XpetraCrsGraphAdapter
52 #include <string>
53 
56 #include <Zoltan2_TestHelpers.hpp>
57 
58 #include <Teuchos_DefaultComm.hpp>
59 #include <Teuchos_RCP.hpp>
60 #include <Teuchos_Comm.hpp>
61 #include <Teuchos_CommHelpers.hpp>
62 
63 using Teuchos::RCP;
64 using Teuchos::rcp;
65 using Teuchos::rcp_const_cast;
66 using Teuchos::Comm;
67 using Teuchos::Array;
68 using Teuchos::ArrayView;
69 
70 typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t> tgraph_t;
71 typedef Xpetra::CrsGraph<zlno_t, zgno_t, znode_t> xgraph_t;
72 
73 template<typename offset_t>
74 void printGraph(RCP<const Comm<int> > &comm, zlno_t nvtx,
75  const zgno_t *vtxIds, const offset_t *offsets, const zgno_t *edgeIds)
76 {
77  int rank = comm->getRank();
78  int nprocs = comm->getSize();
79  comm->barrier();
80  for (int p=0; p < nprocs; p++){
81  if (p == rank){
82  std::cout << rank << ":" << std::endl;
83  for (zlno_t i=0; i < nvtx; i++){
84  std::cout << " vertex " << vtxIds[i] << ": ";
85  for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
86  std::cout << edgeIds[j] << " ";
87  }
88  std::cout << std::endl;
89  }
90  std::cout.flush();
91  }
92  comm->barrier();
93  }
94  comm->barrier();
95 }
96 
97 template <typename User>
100  tgraph_t &graph
101 )
102 {
103  typedef typename Zoltan2::InputTraits<User>::offset_t offset_t;
104 
105  RCP<const Comm<int> > comm = graph.getComm();
106  int fail = 0, gfail=0;
107 
108  if (!fail &&
109  ia.getLocalNumVertices() != graph.getNodeNumRows())
110  fail = 4;
111 
112  if (!fail &&
113  ia.getLocalNumEdges() != graph.getNodeNumEntries())
114  fail = 6;
115 
116  gfail = globalFail(comm, fail);
117 
118  const zgno_t *vtxIds=NULL, *edgeIds=NULL;
119  const offset_t *offsets=NULL;
120  size_t nvtx=0;
121 
122  if (!gfail){
123 
124  nvtx = ia.getLocalNumVertices();
125  ia.getVertexIDsView(vtxIds);
126  ia.getEdgesView(offsets, edgeIds);
127 
128  if (nvtx != graph.getNodeNumRows())
129  fail = 8;
130 
131  gfail = globalFail(comm, fail);
132 
133  if (gfail == 0){
134  printGraph<offset_t>(comm, nvtx, vtxIds, offsets, edgeIds);
135  }
136  else{
137  if (!fail) fail = 10;
138  }
139  }
140  return fail;
141 }
142 
143 int main(int narg, char *arg[])
144 {
145  Tpetra::ScopeGuard tscope(&narg, &arg);
146  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
147 
148  int rank = comm->getRank();
149  int fail = 0, gfail=0;
150  bool aok = true;
151 
152  // Create an object that can give us test Tpetra, Xpetra
153  // and Epetra graphs for testing.
154 
155  RCP<UserInputForTests> uinput;
156 
157  try{
158  uinput =
159  rcp(new UserInputForTests(testDataFilePath,std::string("simple"), comm, true));
160  }
161  catch(std::exception &e){
162  aok = false;
163  std::cout << e.what() << std::endl;
164  }
165  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
166 
167  RCP<tgraph_t> tG; // original graph (for checking)
168  RCP<tgraph_t> newG; // migrated graph
169 
170  tG = uinput->getUITpetraCrsGraph();
171  size_t nvtx = tG->getNodeNumRows();
172 
173  // To test migration in the input adapter we need a Solution object.
174  // Our solution just assigns all objects to part zero.
175 
176  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
177 
178  int nWeights = 1;
179 
182  typedef adapter_t::part_t part_t;
183 
184  part_t *p = new part_t [nvtx];
185  memset(p, 0, sizeof(part_t) * nvtx);
186  ArrayRCP<part_t> solnParts(p, 0, nvtx, true);
187 
188  soln_t solution(env, comm, nWeights);
189  solution.setParts(solnParts);
190 
192  // User object is Tpetra::CrsGraph
193  if (!gfail){
194  if (rank==0)
195  std::cout << "Input adapter for Tpetra::CrsGraph" << std::endl;
196 
197  RCP<const tgraph_t> ctG = rcp_const_cast<const tgraph_t>(tG);
198  RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > tGInput;
199 
200  try {
201  tGInput =
203  }
204  catch (std::exception &e){
205  aok = false;
206  std::cout << e.what() << std::endl;
207  }
208  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter ", 1);
209 
210  fail = verifyInputAdapter<tgraph_t>(*tGInput, *tG);
211 
212  gfail = globalFail(comm, fail);
213 
214  if (!gfail){
215  tgraph_t *mMigrate = NULL;
216  try{
217  tGInput->applyPartitioningSolution( *tG, mMigrate, solution);
218  newG = rcp(mMigrate);
219  }
220  catch (std::exception &e){
221  fail = 11;
222  }
223 
224  gfail = globalFail(comm, fail);
225 
226  if (!gfail){
227  RCP<const tgraph_t> cnewG = rcp_const_cast<const tgraph_t>(newG);
228  RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > newInput;
229  try{
230  newInput = rcp(new Zoltan2::XpetraCrsGraphAdapter<tgraph_t>(cnewG));
231  }
232  catch (std::exception &e){
233  aok = false;
234  std::cout << e.what() << std::endl;
235  }
236  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 2 ", 1);
237 
238  if (rank==0){
239  std::cout <<
240  "Input adapter for Tpetra::CrsGraph migrated to proc 0" <<
241  std::endl;
242  }
243  fail = verifyInputAdapter<tgraph_t>(*newInput, *newG);
244  if (fail) fail += 100;
245  gfail = globalFail(comm, fail);
246  }
247  }
248  if (gfail){
249  printFailureCode(comm, fail);
250  }
251  }
252 
254  // User object is Xpetra::CrsGraph
255  if (!gfail){
256  if (rank==0)
257  std::cout << "Input adapter for Xpetra::CrsGraph" << std::endl;
258 
259  RCP<xgraph_t> xG = uinput->getUIXpetraCrsGraph();
260  RCP<const xgraph_t> cxG = rcp_const_cast<const xgraph_t>(xG);
261  RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > xGInput;
262 
263  try {
264  xGInput =
266  }
267  catch (std::exception &e){
268  aok = false;
269  std::cout << e.what() << std::endl;
270  }
271  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 3 ", 1);
272 
273  fail = verifyInputAdapter<xgraph_t>(*xGInput, *tG);
274 
275  gfail = globalFail(comm, fail);
276 
277  if (!gfail){
278  xgraph_t *mMigrate =NULL;
279  try{
280  xGInput->applyPartitioningSolution( *xG, mMigrate, solution);
281  }
282  catch (std::exception &e){
283  fail = 11;
284  }
285 
286  gfail = globalFail(comm, fail);
287 
288  if (!gfail){
289  RCP<const xgraph_t> cnewG(mMigrate);
290  RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > newInput;
291  try{
292  newInput =
294  }
295  catch (std::exception &e){
296  aok = false;
297  std::cout << e.what() << std::endl;
298  }
299  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 4 ", 1);
300 
301  if (rank==0){
302  std::cout <<
303  "Input adapter for Xpetra::CrsGraph migrated to proc 0" <<
304  std::endl;
305  }
306  fail = verifyInputAdapter<xgraph_t>(*newInput, *newG);
307  if (fail) fail += 100;
308  gfail = globalFail(comm, fail);
309  }
310  }
311  if (gfail){
312  printFailureCode(comm, fail);
313  }
314  }
315 
316 #ifdef HAVE_EPETRA_DATA_TYPES
317  // User object is Epetra_CrsGraph
319  typedef Epetra_CrsGraph egraph_t;
320  if (!gfail){
321  if (rank==0)
322  std::cout << "Input adapter for Epetra_CrsGraph" << std::endl;
323 
324  RCP<egraph_t> eG = uinput->getUIEpetraCrsGraph();
325  RCP<const egraph_t> ceG = rcp_const_cast<const egraph_t>(eG);
326  RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > eGInput;
327 
328  bool goodAdapter = true;
329  try {
330  eGInput =
332  }
333  catch (std::exception &e){
334  if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
335  aok = false;
336  goodAdapter = false;
337  std::cout << e.what() << std::endl;
338  }
339  else {
340  // We expect an error from Xpetra when znode_t != Xpetra::EpetraNode
341  // Ignore it, but skip tests using this matrix adapter.
342  std::cout << "Node type is not supported by Xpetra's Epetra interface;"
343  << " Skipping this test." << std::endl;
344  std::cout << "FYI: Here's the exception message: " << std::endl
345  << e.what() << std::endl;
346  goodAdapter = false;
347  }
348  }
349  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraCrsGraphAdapter 5 ", 1);
350 
351  if (goodAdapter) {
352  fail = verifyInputAdapter<egraph_t>(*eGInput, *tG);
353 
354  gfail = globalFail(comm, fail);
355 
356  if (!gfail){
357  egraph_t *mMigrate =NULL;
358  try{
359  eGInput->applyPartitioningSolution( *eG, mMigrate, solution);
360  }
361  catch (std::exception &e){
362  fail = 11;
363  }
364 
365  gfail = globalFail(comm, fail);
366 
367  if (!gfail){
368  RCP<const egraph_t> cnewG(mMigrate, true);
369  RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_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, "XpetraCrsGraphAdapter 6 ", 1);
379 
380  if (rank==0){
381  std::cout <<
382  "Input adapter for Epetra_CrsGraph migrated to proc 0" <<
383  std::endl;
384  }
385  fail = verifyInputAdapter<egraph_t>(*newInput, *newG);
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::XpetraCrsGraphAdapter
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Definition: Zoltan2_XpetraCrsGraphAdapter.hpp:84
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
printGraph
void printGraph(RCP< const Comm< int > > &comm, zlno_t nvtx, const zgno_t *vtxIds, const offset_t *offsets, const zgno_t *edgeIds)
Definition: XpetraCrsGraphInput.cpp:74
Zoltan2::InputTraits::offset_t
default_offset_t offset_t
The data type to represent offsets.
Definition: Zoltan2_InputTraits.hpp:195
testDataFilePath
std::string testDataFilePath(".")
UserInputForTests
Definition: UserInputForTests.hpp:126
part_t
SparseMatrixAdapter_t::part_t part_t
Definition: partitioningTree.cpp:74
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::XpetraCrsGraphAdapter::getLocalNumVertices
size_t getLocalNumVertices() const
Returns the number of vertices on this process.
Definition: Zoltan2_XpetraCrsGraphAdapter.hpp:203
Zoltan2::XpetraCrsGraphAdapter::getEdgesView
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
Definition: Zoltan2_XpetraCrsGraphAdapter.hpp:214
Zoltan2_PartitioningSolution.hpp
Defines the PartitioningSolution class.
verifyInputAdapter
int verifyInputAdapter(Zoltan2::XpetraCrsGraphAdapter< User > &ia, tgraph_t &graph)
Definition: XpetraCrsGraphInput.cpp:98
Zoltan2::XpetraCrsGraphAdapter::getLocalNumEdges
size_t getLocalNumEdges() const
Returns the number of edges on this process.
Definition: Zoltan2_XpetraCrsGraphAdapter.hpp:212
zlno_t
int zlno_t
Definition: Zoltan2_TestHelpers.hpp:142
tgraph_t
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tgraph_t
Definition: XpetraCrsGraphInput.cpp:70
Zoltan2_XpetraCrsGraphAdapter.hpp
Defines XpetraCrsGraphAdapter class.
fail
static const std::string fail
Definition: findUniqueGids.cpp:80
Zoltan2_TestHelpers.hpp
common code used by tests
Zoltan2::XpetraCrsGraphAdapter::getVertexIDsView
void getVertexIDsView(const gno_t *&ids) const
Sets pointers to this process' graph entries.
Definition: Zoltan2_XpetraCrsGraphAdapter.hpp:205
xgraph_t
Xpetra::CrsGraph< zlno_t, zgno_t, znode_t > xgraph_t
Definition: XpetraCrsGraphInput.cpp:71
main
int main(int narg, char *arg[])
Definition: XpetraCrsGraphInput.cpp:143