Zoltan2
TpetraRowGraphInput.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::TpetraRowGraphAdapter
52 #include <string>
53 
55 #include <Zoltan2_InputTraits.hpp>
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::rcp_dynamic_cast;
67 using Teuchos::Comm;
68 
69 typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t> ztcrsgraph_t;
70 typedef Tpetra::RowGraph<zlno_t, zgno_t, znode_t> ztrowgraph_t;
71 
72 template<typename offset_t>
73 void printGraph(RCP<const Comm<int> > &comm, zlno_t nvtx,
74  const zgno_t *vtxIds, const offset_t *offsets, const zgno_t *edgeIds)
75 {
76  int rank = comm->getRank();
77  int nprocs = comm->getSize();
78  comm->barrier();
79  for (int p=0; p < nprocs; p++){
80  if (p == rank){
81  std::cout << rank << ":" << std::endl;
82  for (zlno_t i=0; i < nvtx; i++){
83  std::cout << " vertex " << vtxIds[i] << ": ";
84  for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
85  std::cout << edgeIds[j] << " ";
86  }
87  std::cout << std::endl;
88  }
89  std::cout.flush();
90  }
91  comm->barrier();
92  }
93  comm->barrier();
94 }
95 
96 template <typename User>
99 {
100  typedef typename Zoltan2::InputTraits<User>::offset_t offset_t;
101 
102  RCP<const Comm<int> > comm = graph.getComm();
103  int fail = 0, gfail=0;
104 
105  if (!fail &&
106  ia.getLocalNumVertices() != graph.getNodeNumRows())
107  fail = 4;
108 
109  if (!fail &&
110  ia.getLocalNumEdges() != graph.getNodeNumEntries())
111  fail = 6;
112 
113  gfail = globalFail(comm, fail);
114 
115  const zgno_t *vtxIds=NULL, *edgeIds=NULL;
116  const offset_t *offsets=NULL;
117  size_t nvtx=0;
118 
119  if (!gfail){
120 
121  nvtx = ia.getLocalNumVertices();
122  ia.getVertexIDsView(vtxIds);
123  ia.getEdgesView(offsets, edgeIds);
124 
125  if (nvtx != graph.getNodeNumRows())
126  fail = 8;
127 
128  gfail = globalFail(comm, fail);
129 
130  if (gfail == 0){
131  printGraph<offset_t>(comm, nvtx, vtxIds, offsets, edgeIds);
132  }
133  else{
134  if (!fail) fail = 10;
135  }
136  }
137  return fail;
138 }
139 
140 int main(int narg, char *arg[])
141 {
142  Tpetra::ScopeGuard tscope(&narg, &arg);
143  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
144 
145  int rank = comm->getRank();
146  int fail = 0, gfail=0;
147  bool aok = true;
148 
149  // Create an object that can give us test Tpetra graphs for testing
150 
151  RCP<UserInputForTests> uinput;
152 
153  try{
154  uinput =
155  rcp(new UserInputForTests(
156  testDataFilePath,std::string("simple"), comm, true));
157  }
158  catch(std::exception &e){
159  aok = false;
160  std::cout << e.what() << std::endl;
161  }
162  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
163 
164  // Input crs graph and row graph cast from it.
165  RCP<ztcrsgraph_t> tG = uinput->getUITpetraCrsGraph();
166  RCP<ztrowgraph_t> trG = rcp_dynamic_cast<ztrowgraph_t>(tG);
167 
168  RCP<ztrowgraph_t> newG; // migrated graph
169 
170  size_t nvtx = tG->getNodeNumRows();
171 
172  // To test migration in the input adapter we need a Solution object.
173 
174  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
175 
176  int nWeights = 1;
177 
180  typedef adapter_t::part_t part_t;
181 
182  part_t *p = new part_t [nvtx];
183  memset(p, 0, sizeof(part_t) * nvtx);
184  ArrayRCP<part_t> solnParts(p, 0, nvtx, true);
185 
186  soln_t solution(env, comm, nWeights);
187  solution.setParts(solnParts);
188 
190  // User object is Tpetra::RowGraph
191  if (!gfail){
192  if (rank==0)
193  std::cout << "Input adapter for Tpetra::RowGraph" << std::endl;
194 
195  RCP<const ztrowgraph_t> ctrG = rcp_const_cast<const ztrowgraph_t>(
196  rcp_dynamic_cast<ztrowgraph_t>(tG));
197 
198  RCP<adapter_t> trGInput;
199 
200  try {
201  trGInput = rcp(new adapter_t(ctrG));
202  }
203  catch (std::exception &e){
204  aok = false;
205  std::cout << e.what() << std::endl;
206  }
207  TEST_FAIL_AND_EXIT(*comm, aok, "TpetraRowGraphAdapter ", 1);
208 
209  fail = verifyInputAdapter<ztrowgraph_t>(*trGInput, *trG);
210 
211  gfail = globalFail(comm, fail);
212 
213  if (!gfail){
214  ztrowgraph_t *mMigrate = NULL;
215  try{
216  trGInput->applyPartitioningSolution( *trG, mMigrate, solution);
217  newG = rcp(mMigrate);
218  }
219  catch (std::exception &e){
220  fail = 11;
221  }
222 
223  gfail = globalFail(comm, fail);
224 
225  if (!gfail){
226  RCP<const ztrowgraph_t> cnewG =
227  rcp_const_cast<const ztrowgraph_t>(newG);
228  RCP<adapter_t> newInput;
229  try{
230  newInput = rcp(new adapter_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, "TpetraRowGraphAdapter 2 ", 1);
237 
238  if (rank==0){
239  std::cout <<
240  "Input adapter for Tpetra::RowGraph migrated to proc 0" <<
241  std::endl;
242  }
243  fail = verifyInputAdapter<ztrowgraph_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  // DONE
255 
256  if (rank==0)
257  std::cout << "PASS" << std::endl;
258 }
Zoltan2_InputTraits.hpp
Traits for application input objects.
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
Zoltan2_TpetraRowGraphAdapter.hpp
Defines TpetraRowGraphAdapter class.
Zoltan2::InputTraits::offset_t
default_offset_t offset_t
The data type to represent offsets.
Definition: Zoltan2_InputTraits.hpp:195
testDataFilePath
std::string testDataFilePath(".")
Zoltan2::TpetraRowGraphAdapter
Provides access for Zoltan2 to Tpetra::RowGraph data.
Definition: Zoltan2_TpetraRowGraphAdapter.hpp:82
printGraph
void printGraph(RCP< const Comm< int > > &comm, zlno_t nvtx, const zgno_t *vtxIds, const offset_t *offsets, const zgno_t *edgeIds)
Definition: TpetraRowGraphInput.cpp:73
ztcrsgraph_t
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > ztcrsgraph_t
Definition: TpetraRowGraphInput.cpp:69
UserInputForTests
Definition: UserInputForTests.hpp:126
part_t
SparseMatrixAdapter_t::part_t part_t
Definition: partitioningTree.cpp:74
Zoltan2::TpetraRowGraphAdapter::getVertexIDsView
void getVertexIDsView(const gno_t *&ids) const
Sets pointers to this process' graph entries.
Definition: Zoltan2_TpetraRowGraphAdapter.hpp:194
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
ztrowgraph_t
Tpetra::RowGraph< zlno_t, zgno_t, znode_t > ztrowgraph_t
Definition: TpetraRowGraphInput.cpp:70
Zoltan2::TpetraRowGraphAdapter::getLocalNumEdges
size_t getLocalNumEdges() const
Returns the number of edges on this process.
Definition: Zoltan2_TpetraRowGraphAdapter.hpp:201
Zoltan2::TpetraRowGraphAdapter::getEdgesView
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const
Definition: Zoltan2_TpetraRowGraphAdapter.hpp:203
zlno_t
int zlno_t
Definition: Zoltan2_TestHelpers.hpp:142
main
int main(int narg, char *arg[])
Definition: TpetraRowGraphInput.cpp:140
verifyInputAdapter
int verifyInputAdapter(Zoltan2::TpetraRowGraphAdapter< User > &ia, ztrowgraph_t &graph)
Definition: TpetraRowGraphInput.cpp:97
fail
static const std::string fail
Definition: findUniqueGids.cpp:80
Zoltan2_TestHelpers.hpp
common code used by tests
Zoltan2::TpetraRowGraphAdapter::getLocalNumVertices
size_t getLocalNumVertices() const
Returns the number of vertices on this process.
Definition: Zoltan2_TpetraRowGraphAdapter.hpp:192