Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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  auto comm = graph.getComm();
103  int fail = 0, gfail=0;
104 
105  if (!fail &&
106  ia.getLocalNumVertices() != graph.getLocalNumRows())
107  fail = 4;
108 
109  if (!fail &&
110  ia.getLocalNumEdges() != graph.getLocalNumEntries())
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.getLocalNumRows())
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  Teuchos::ParameterList params;
153  params.set("input file", "simple");
154  params.set("file type", "Chaco");
155 
156  try{
157  uinput = rcp(new UserInputForTests(params, comm));
158  }
159  catch(std::exception &e){
160  aok = false;
161  std::cout << e.what() << std::endl;
162  }
163  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
164 
165  // Input crs graph and row graph cast from it.
166  auto tG = uinput->getUITpetraCrsGraph();
167  auto trG = rcp_dynamic_cast<ztrowgraph_t>(tG);
168 
169  RCP<ztrowgraph_t> newG; // migrated graph
170 
171  size_t nvtx = tG->getLocalNumRows();
172 
173  // To test migration in the input adapter we need a Solution object.
174 
175  const auto env = rcp(new Zoltan2::Environment(comm));
176 
177  int nWeights = 1;
178 
181  typedef adapter_t::part_t part_t;
182 
183  part_t *p = new part_t [nvtx];
184  memset(p, 0, sizeof(part_t) * nvtx);
185  ArrayRCP<part_t> solnParts(p, 0, nvtx, true);
186 
187  soln_t solution(env, comm, nWeights);
188  solution.setParts(solnParts);
189 
191  // User object is Tpetra::RowGraph
192  if (!gfail){
193  if (rank==0)
194  std::cout << "Input adapter for Tpetra::RowGraph" << std::endl;
195 
196  RCP<const ztrowgraph_t> ctrG = rcp_const_cast<const ztrowgraph_t>(
197  rcp_dynamic_cast<ztrowgraph_t>(tG));
198 
199  RCP<adapter_t> trGInput;
200 
201  try {
202  trGInput = rcp(new adapter_t(ctrG));
203  }
204  catch (std::exception &e){
205  aok = false;
206  std::cout << e.what() << std::endl;
207  }
208  TEST_FAIL_AND_EXIT(*comm, aok, "TpetraRowGraphAdapter ", 1);
209 
210  fail = verifyInputAdapter<ztrowgraph_t>(*trGInput, *trG);
211 
212  gfail = globalFail(*comm, fail);
213 
214  if (!gfail){
215  ztrowgraph_t *mMigrate = NULL;
216  try{
217  trGInput->applyPartitioningSolution( *trG, 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  auto cnewG = 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 }
void printFailureCode(const Comm< int > &comm, int fail)
void verifyInputAdapter(adapter_t &ia, matrix_t &matrix)
void getVertexIDsView(const gno_t *&ids) const override
Sets pointers to this process&#39; graph entries.
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const override
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
Defines TpetraRowGraphAdapter class.
default_offset_t offset_t
The data type to represent offsets.
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
size_t getLocalNumVertices() const override
Returns the number of vertices on this process.
int main(int narg, char **arg)
Definition: coloring1.cpp:199
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)
SparseMatrixAdapter_t::part_t part_t
Tpetra::RowGraph< zlno_t, zgno_t, znode_t > ztrowgraph_t
A PartitioningSolution is a solution to a partitioning problem.
size_t getLocalNumEdges() const override
Returns the number of edges on this process.
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > ztcrsgraph_t
Traits for application input objects.
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)
Tpetra::Map::global_ordinal_type zgno_t
Provides access for Zoltan2 to Tpetra::RowGraph data.