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