Zoltan2
Loading...
Searching...
No Matches
TpetraRowGraphInput.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::TpetraRowGraphAdapter
17#include <string>
18
22
23#include <Teuchos_DefaultComm.hpp>
24#include <Teuchos_RCP.hpp>
25#include <Teuchos_Comm.hpp>
26#include <Teuchos_CommHelpers.hpp>
27
28using Teuchos::RCP;
29using Teuchos::rcp;
30using Teuchos::rcp_const_cast;
31using Teuchos::rcp_dynamic_cast;
32using Teuchos::Comm;
33
34typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t> ztcrsgraph_t;
35typedef Tpetra::RowGraph<zlno_t, zgno_t, znode_t> ztrowgraph_t;
36
37template<typename offset_t>
38void printGraph(RCP<const Comm<int> > &comm, zlno_t nvtx,
39 const zgno_t *vtxIds, const offset_t *offsets, const zgno_t *edgeIds)
40{
41 int rank = comm->getRank();
42 int nprocs = comm->getSize();
43 comm->barrier();
44 for (int p=0; p < nprocs; p++){
45 if (p == rank){
46 std::cout << rank << ":" << std::endl;
47 for (zlno_t i=0; i < nvtx; i++){
48 std::cout << " vertex " << vtxIds[i] << ": ";
49 for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
50 std::cout << edgeIds[j] << " ";
51 }
52 std::cout << std::endl;
53 }
54 std::cout.flush();
55 }
56 comm->barrier();
57 }
58 comm->barrier();
59}
60
61template <typename User>
64{
65 typedef typename Zoltan2::InputTraits<User>::offset_t offset_t;
66
67 auto comm = graph.getComm();
68 int fail = 0, gfail=0;
69
70 if (!fail &&
71 ia.getLocalNumVertices() != graph.getLocalNumRows())
72 fail = 4;
73
74 if (!fail &&
75 ia.getLocalNumEdges() != graph.getLocalNumEntries())
76 fail = 6;
77
78 gfail = globalFail(*comm, fail);
79
80 const zgno_t *vtxIds=NULL, *edgeIds=NULL;
81 const offset_t *offsets=NULL;
82 size_t nvtx=0;
83
84 if (!gfail){
85
86 nvtx = ia.getLocalNumVertices();
87 ia.getVertexIDsView(vtxIds);
88 ia.getEdgesView(offsets, edgeIds);
89
90 if (nvtx != graph.getLocalNumRows())
91 fail = 8;
92
93 gfail = globalFail(*comm, fail);
94
95 if (gfail == 0){
96 printGraph<offset_t>(comm, nvtx, vtxIds, offsets, edgeIds);
97 }
98 else{
99 if (!fail) fail = 10;
100 }
101 }
102 return fail;
103}
104
105int main(int narg, char *arg[])
106{
107 Tpetra::ScopeGuard tscope(&narg, &arg);
108 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
109
110 int rank = comm->getRank();
111 int fail = 0, gfail=0;
112 bool aok = true;
113
114 // Create an object that can give us test Tpetra graphs for testing
115
116 RCP<UserInputForTests> uinput;
117 Teuchos::ParameterList params;
118 params.set("input file", "simple");
119 params.set("file type", "Chaco");
120
121 try{
122 uinput = rcp(new UserInputForTests(params, comm));
123 }
124 catch(std::exception &e){
125 aok = false;
126 std::cout << e.what() << std::endl;
127 }
128 TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
129
130 // Input crs graph and row graph cast from it.
131 auto tG = uinput->getUITpetraCrsGraph();
132 auto trG = rcp_dynamic_cast<ztrowgraph_t>(tG);
133
134 RCP<ztrowgraph_t> newG; // migrated graph
135
136 size_t nvtx = tG->getLocalNumRows();
137
138 // To test migration in the input adapter we need a Solution object.
139
140 const auto env = rcp(new Zoltan2::Environment(comm));
141
142 int nWeights = 1;
143
146 typedef adapter_t::part_t part_t;
147
148 part_t *p = new part_t [nvtx];
149 memset(p, 0, sizeof(part_t) * nvtx);
150 ArrayRCP<part_t> solnParts(p, 0, nvtx, true);
151
152 soln_t solution(env, comm, nWeights);
153 solution.setParts(solnParts);
154
156 // User object is Tpetra::RowGraph
157 if (!gfail){
158 if (rank==0)
159 std::cout << "Input adapter for Tpetra::RowGraph" << std::endl;
160
161 RCP<const ztrowgraph_t> ctrG = rcp_const_cast<const ztrowgraph_t>(
162 rcp_dynamic_cast<ztrowgraph_t>(tG));
163
164 RCP<adapter_t> trGInput;
165
166 try {
167 trGInput = rcp(new adapter_t(ctrG));
168 }
169 catch (std::exception &e){
170 aok = false;
171 std::cout << e.what() << std::endl;
172 }
173 TEST_FAIL_AND_EXIT(*comm, aok, "TpetraRowGraphAdapter ", 1);
174
175 fail = verifyInputAdapter<ztrowgraph_t>(*trGInput, *trG);
176
177 gfail = globalFail(*comm, fail);
178
179 if (!gfail){
180 ztrowgraph_t *mMigrate = NULL;
181 try{
182 trGInput->applyPartitioningSolution( *trG, mMigrate, solution);
183 newG = rcp(mMigrate);
184 }
185 catch (std::exception &e){
186 fail = 11;
187 }
188
189 gfail = globalFail(*comm, fail);
190
191 if (!gfail){
192 auto cnewG = rcp_const_cast<const ztrowgraph_t>(newG);
193 RCP<adapter_t> newInput;
194 try{
195 newInput = rcp(new adapter_t(cnewG));
196 }
197 catch (std::exception &e){
198 aok = false;
199 std::cout << e.what() << std::endl;
200 }
201 TEST_FAIL_AND_EXIT(*comm, aok, "TpetraRowGraphAdapter 2 ", 1);
202
203 if (rank==0){
204 std::cout <<
205 "Input adapter for Tpetra::RowGraph migrated to proc 0" <<
206 std::endl;
207 }
208 fail = verifyInputAdapter<ztrowgraph_t>(*newInput, *newG);
209 if (fail) fail += 100;
210 gfail = globalFail(*comm, fail);
211 }
212 }
213 if (gfail){
214 printFailureCode(*comm, fail);
215 }
216 }
217
219 // DONE
220
221 if (rank==0)
222 std::cout << "PASS" << std::endl;
223}
int globalFail(const Comm< int > &comm, int fail)
void printFailureCode(const Comm< int > &comm, int fail)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Tpetra::RowGraph< zlno_t, zgno_t, znode_t > ztrowgraph_t
void printGraph(RCP< const Comm< int > > &comm, zlno_t nvtx, const zgno_t *vtxIds, const offset_t *offsets, const zgno_t *edgeIds)
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > ztcrsgraph_t
int verifyInputAdapter(Zoltan2::TpetraRowGraphAdapter< User > &ia, ztrowgraph_t &graph)
Traits for application input objects.
common code used by tests
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
Defines TpetraRowGraphAdapter class.
int main()
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
A PartitioningSolution is a solution to a partitioning problem.
Provides access for Zoltan2 to Tpetra::RowGraph data.
size_t getLocalNumVertices() const override
Returns the number of vertices on this process.
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.
void getVertexIDsView(const gno_t *&ids) const override
Sets pointers to this process' graph entries.
size_t getLocalNumEdges() const override
Returns the number of edges on this process.
static const std::string fail
SparseMatrixAdapter_t::part_t part_t
default_offset_t offset_t
The data type to represent offsets.