58 #include <Teuchos_DefaultComm.hpp>
59 #include <Teuchos_RCP.hpp>
60 #include <Teuchos_Comm.hpp>
61 #include <Teuchos_CommHelpers.hpp>
65 using Teuchos::rcp_const_cast;
68 using Teuchos::ArrayView;
70 typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t>
tgraph_t;
71 typedef Xpetra::CrsGraph<zlno_t, zgno_t, znode_t>
xgraph_t;
73 template<
typename offset_t>
75 const zgno_t *vtxIds,
const offset_t *offsets,
const zgno_t *edgeIds)
77 int rank = comm->getRank();
78 int nprocs = comm->getSize();
80 for (
int p=0; p < nprocs; p++){
82 std::cout << rank <<
":" << std::endl;
83 for (
zlno_t i=0; i < nvtx; i++){
84 std::cout <<
" vertex " << vtxIds[i] <<
": ";
85 for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
86 std::cout << edgeIds[j] <<
" ";
88 std::cout << std::endl;
97 template <
typename User>
105 RCP<const Comm<int> > comm = graph.getComm();
106 int fail = 0, gfail=0;
118 const zgno_t *vtxIds=NULL, *edgeIds=NULL;
119 const offset_t *offsets=NULL;
128 if (nvtx != graph.getNodeNumRows())
134 printGraph<offset_t>(comm, nvtx, vtxIds, offsets, edgeIds);
143 int main(
int narg,
char *arg[])
145 Tpetra::ScopeGuard tscope(&narg, &arg);
146 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
148 int rank = comm->getRank();
149 int fail = 0, gfail=0;
155 RCP<UserInputForTests> uinput;
161 catch(std::exception &e){
163 std::cout << e.what() << std::endl;
170 tG = uinput->getUITpetraCrsGraph();
171 size_t nvtx = tG->getNodeNumRows();
185 memset(p, 0,
sizeof(
part_t) * nvtx);
186 ArrayRCP<part_t> solnParts(p, 0, nvtx,
true);
188 soln_t solution(env, comm, nWeights);
189 solution.setParts(solnParts);
195 std::cout <<
"Input adapter for Tpetra::CrsGraph" << std::endl;
197 RCP<const tgraph_t> ctG = rcp_const_cast<const tgraph_t>(tG);
198 RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > tGInput;
204 catch (std::exception &e){
206 std::cout << e.what() << std::endl;
210 fail = verifyInputAdapter<tgraph_t>(*tGInput, *tG);
217 tGInput->applyPartitioningSolution( *tG, mMigrate, solution);
218 newG = rcp(mMigrate);
220 catch (std::exception &e){
227 RCP<const tgraph_t> cnewG = rcp_const_cast<const tgraph_t>(newG);
228 RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > newInput;
232 catch (std::exception &e){
234 std::cout << e.what() << std::endl;
240 "Input adapter for Tpetra::CrsGraph migrated to proc 0" <<
243 fail = verifyInputAdapter<tgraph_t>(*newInput, *newG);
257 std::cout <<
"Input adapter for Xpetra::CrsGraph" << std::endl;
259 RCP<xgraph_t> xG = uinput->getUIXpetraCrsGraph();
260 RCP<const xgraph_t> cxG = rcp_const_cast<const xgraph_t>(xG);
261 RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > xGInput;
267 catch (std::exception &e){
269 std::cout << e.what() << std::endl;
273 fail = verifyInputAdapter<xgraph_t>(*xGInput, *tG);
280 xGInput->applyPartitioningSolution( *xG, mMigrate, solution);
282 catch (std::exception &e){
289 RCP<const xgraph_t> cnewG(mMigrate);
290 RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > newInput;
295 catch (std::exception &e){
297 std::cout << e.what() << std::endl;
303 "Input adapter for Xpetra::CrsGraph migrated to proc 0" <<
306 fail = verifyInputAdapter<xgraph_t>(*newInput, *newG);
316 #ifdef HAVE_EPETRA_DATA_TYPES
319 typedef Epetra_CrsGraph egraph_t;
322 std::cout <<
"Input adapter for Epetra_CrsGraph" << std::endl;
324 RCP<egraph_t> eG = uinput->getUIEpetraCrsGraph();
325 RCP<const egraph_t> ceG = rcp_const_cast<const egraph_t>(eG);
326 RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > eGInput;
328 bool goodAdapter =
true;
333 catch (std::exception &e){
334 if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
337 std::cout << e.what() << std::endl;
342 std::cout <<
"Node type is not supported by Xpetra's Epetra interface;"
343 <<
" Skipping this test." << std::endl;
344 std::cout <<
"FYI: Here's the exception message: " << std::endl
345 << e.what() << std::endl;
352 fail = verifyInputAdapter<egraph_t>(*eGInput, *tG);
357 egraph_t *mMigrate =NULL;
359 eGInput->applyPartitioningSolution( *eG, mMigrate, solution);
361 catch (std::exception &e){
368 RCP<const egraph_t> cnewG(mMigrate,
true);
369 RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > newInput;
374 catch (std::exception &e){
376 std::cout << e.what() << std::endl;
382 "Input adapter for Epetra_CrsGraph migrated to proc 0" <<
385 fail = verifyInputAdapter<egraph_t>(*newInput, *newG);
401 std::cout <<
"PASS" << std::endl;