23 #include <Teuchos_DefaultComm.hpp>
24 #include <Teuchos_RCP.hpp>
25 #include <Teuchos_Comm.hpp>
26 #include <Teuchos_CommHelpers.hpp>
30 using Teuchos::rcp_const_cast;
33 using Teuchos::ArrayView;
35 typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t>
tgraph_t;
36 typedef Xpetra::CrsGraph<zlno_t, zgno_t, znode_t>
xgraph_t;
38 template<
typename offset_t>
40 const zgno_t *vtxIds,
const offset_t *offsets,
const zgno_t *edgeIds)
42 int rank = comm->getRank();
43 int nprocs = comm->getSize();
45 for (
int p=0; p < nprocs; p++){
47 std::cout << rank <<
":" << std::endl;
48 for (
zlno_t i=0; i < nvtx; i++){
49 std::cout <<
" vertex " << vtxIds[i] <<
": ";
50 for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
51 std::cout << edgeIds[j] <<
" ";
53 std::cout << std::endl;
62 template <
typename User>
70 RCP<const Comm<int> > comm = graph.getComm();
71 int fail = 0, gfail=0;
83 const zgno_t *vtxIds=NULL, *edgeIds=NULL;
84 const offset_t *offsets=NULL;
93 if (nvtx != graph.getLocalNumRows())
99 printGraph<offset_t>(comm, nvtx, vtxIds, offsets, edgeIds);
102 if (!fail) fail = 10;
108 int main(
int narg,
char *arg[])
110 Tpetra::ScopeGuard tscope(&narg, &arg);
111 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
113 int rank = comm->getRank();
114 int fail = 0, gfail=0;
120 RCP<UserInputForTests> uinput;
121 Teuchos::ParameterList params;
122 params.set(
"input file",
"simple");
123 params.set(
"file type",
"Chaco");
128 catch(std::exception &e){
130 std::cout << e.what() << std::endl;
137 tG = uinput->getUITpetraCrsGraph();
138 size_t nvtx = tG->getLocalNumRows();
151 part_t *p =
new part_t [nvtx];
152 memset(p, 0,
sizeof(part_t) * nvtx);
153 ArrayRCP<part_t> solnParts(p, 0, nvtx,
true);
155 soln_t solution(env, comm, nWeights);
156 solution.setParts(solnParts);
162 std::cout <<
"Input adapter for Tpetra::CrsGraph" << std::endl;
164 RCP<const tgraph_t> ctG = rcp_const_cast<
const tgraph_t>(tG);
165 RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > tGInput;
171 catch (std::exception &e){
173 std::cout << e.what() << std::endl;
177 fail = verifyInputAdapter<tgraph_t>(*tGInput, *tG);
184 tGInput->applyPartitioningSolution( *tG, mMigrate, solution);
185 newG = rcp(mMigrate);
187 catch (std::exception &e){
194 RCP<const tgraph_t> cnewG = rcp_const_cast<
const tgraph_t>(newG);
195 RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > newInput;
199 catch (std::exception &e){
201 std::cout << e.what() << std::endl;
207 "Input adapter for Tpetra::CrsGraph migrated to proc 0" <<
210 fail = verifyInputAdapter<tgraph_t>(*newInput, *newG);
211 if (fail) fail += 100;
224 std::cout <<
"Input adapter for Xpetra::CrsGraph" << std::endl;
226 RCP<xgraph_t> xG = uinput->getUIXpetraCrsGraph();
227 RCP<const xgraph_t> cxG = rcp_const_cast<
const xgraph_t>(xG);
228 RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > xGInput;
234 catch (std::exception &e){
236 std::cout << e.what() << std::endl;
240 fail = verifyInputAdapter<xgraph_t>(*xGInput, *tG);
247 xGInput->applyPartitioningSolution( *xG, mMigrate, solution);
249 catch (std::exception &e){
256 RCP<const xgraph_t> cnewG(mMigrate);
257 RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > newInput;
262 catch (std::exception &e){
264 std::cout << e.what() << std::endl;
270 "Input adapter for Xpetra::CrsGraph migrated to proc 0" <<
273 fail = verifyInputAdapter<xgraph_t>(*newInput, *newG);
274 if (fail) fail += 100;
283 #ifdef HAVE_EPETRA_DATA_TYPES
286 typedef Epetra_CrsGraph egraph_t;
289 std::cout <<
"Input adapter for Epetra_CrsGraph" << std::endl;
291 RCP<egraph_t> eG = uinput->getUIEpetraCrsGraph();
292 RCP<const egraph_t> ceG = rcp_const_cast<
const egraph_t>(eG);
293 RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > eGInput;
295 bool goodAdapter =
true;
300 catch (std::exception &e){
301 if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
304 std::cout << e.what() << std::endl;
309 std::cout <<
"Node type is not supported by Xpetra's Epetra interface;"
310 <<
" Skipping this test." << std::endl;
311 std::cout <<
"FYI: Here's the exception message: " << std::endl
312 << e.what() << std::endl;
319 fail = verifyInputAdapter<egraph_t>(*eGInput, *tG);
324 egraph_t *mMigrate =NULL;
326 eGInput->applyPartitioningSolution( *eG, mMigrate, solution);
328 catch (std::exception &e){
335 RCP<const egraph_t> cnewG(mMigrate,
true);
336 RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > newInput;
341 catch (std::exception &e){
343 std::cout << e.what() << std::endl;
349 "Input adapter for Epetra_CrsGraph migrated to proc 0" <<
352 fail = verifyInputAdapter<egraph_t>(*newInput, *newG);
353 if (fail) fail += 100;
368 std::cout <<
"PASS" << std::endl;
void printFailureCode(const Comm< int > &comm, int fail)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Provides access for Zoltan2 to Xpetra::CrsGraph data.
int main(int narg, char **arg)
size_t getLocalNumVertices() const override
Returns the number of vertices on this process.
Defines the PartitioningSolution class.
common code used by tests
Defines XpetraCrsGraphAdapter class.
SparseMatrixAdapter_t::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
size_t getLocalNumEdges() const override
Returns the number of edges on this process.
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)
void getVertexIDsView(const gno_t *&ids) const override
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const override
Tpetra::Map::global_ordinal_type zgno_t