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);
137 if (!fail) fail = 10;
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;
156 Teuchos::ParameterList params;
157 params.set(
"input file",
"simple");
158 params.set(
"file type",
"Chaco");
163 catch(std::exception &e){
165 std::cout << e.what() << std::endl;
172 tG = uinput->getUITpetraCrsGraph();
173 size_t nvtx = tG->getNodeNumRows();
186 part_t *p =
new part_t [nvtx];
187 memset(p, 0,
sizeof(part_t) * nvtx);
188 ArrayRCP<part_t> solnParts(p, 0, nvtx,
true);
190 soln_t solution(env, comm, nWeights);
191 solution.setParts(solnParts);
197 std::cout <<
"Input adapter for Tpetra::CrsGraph" << std::endl;
199 RCP<const tgraph_t> ctG = rcp_const_cast<
const tgraph_t>(tG);
200 RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > tGInput;
206 catch (std::exception &e){
208 std::cout << e.what() << std::endl;
212 fail = verifyInputAdapter<tgraph_t>(*tGInput, *tG);
219 tGInput->applyPartitioningSolution( *tG, mMigrate, solution);
220 newG = rcp(mMigrate);
222 catch (std::exception &e){
229 RCP<const tgraph_t> cnewG = rcp_const_cast<
const tgraph_t>(newG);
230 RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > newInput;
234 catch (std::exception &e){
236 std::cout << e.what() << std::endl;
242 "Input adapter for Tpetra::CrsGraph migrated to proc 0" <<
245 fail = verifyInputAdapter<tgraph_t>(*newInput, *newG);
246 if (fail) fail += 100;
259 std::cout <<
"Input adapter for Xpetra::CrsGraph" << std::endl;
261 RCP<xgraph_t> xG = uinput->getUIXpetraCrsGraph();
262 RCP<const xgraph_t> cxG = rcp_const_cast<
const xgraph_t>(xG);
263 RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > xGInput;
269 catch (std::exception &e){
271 std::cout << e.what() << std::endl;
275 fail = verifyInputAdapter<xgraph_t>(*xGInput, *tG);
282 xGInput->applyPartitioningSolution( *xG, mMigrate, solution);
284 catch (std::exception &e){
291 RCP<const xgraph_t> cnewG(mMigrate);
292 RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > newInput;
297 catch (std::exception &e){
299 std::cout << e.what() << std::endl;
305 "Input adapter for Xpetra::CrsGraph migrated to proc 0" <<
308 fail = verifyInputAdapter<xgraph_t>(*newInput, *newG);
309 if (fail) fail += 100;
318 #ifdef HAVE_EPETRA_DATA_TYPES
321 typedef Epetra_CrsGraph egraph_t;
324 std::cout <<
"Input adapter for Epetra_CrsGraph" << std::endl;
326 RCP<egraph_t> eG = uinput->getUIEpetraCrsGraph();
327 RCP<const egraph_t> ceG = rcp_const_cast<
const egraph_t>(eG);
328 RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > eGInput;
330 bool goodAdapter =
true;
335 catch (std::exception &e){
336 if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
339 std::cout << e.what() << std::endl;
344 std::cout <<
"Node type is not supported by Xpetra's Epetra interface;"
345 <<
" Skipping this test." << std::endl;
346 std::cout <<
"FYI: Here's the exception message: " << std::endl
347 << e.what() << std::endl;
354 fail = verifyInputAdapter<egraph_t>(*eGInput, *tG);
359 egraph_t *mMigrate =NULL;
361 eGInput->applyPartitioningSolution( *eG, mMigrate, solution);
363 catch (std::exception &e){
370 RCP<const egraph_t> cnewG(mMigrate,
true);
371 RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > newInput;
376 catch (std::exception &e){
378 std::cout << e.what() << std::endl;
384 "Input adapter for Epetra_CrsGraph migrated to proc 0" <<
387 fail = verifyInputAdapter<egraph_t>(*newInput, *newG);
388 if (fail) fail += 100;
403 std::cout <<
"PASS" << std::endl;
void printFailureCode(const Comm< int > &comm, int fail)
int main(int narg, char *arg[])
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Defines the PartitioningSolution class.
common code used by tests
size_t getLocalNumEdges() const
Returns the number of edges on this process.
Defines XpetraCrsGraphAdapter class.
SparseMatrixAdapter_t::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
static const std::string fail
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
int globalFail(const Comm< int > &comm, int fail)
void getVertexIDsView(const gno_t *&ids) const
Sets pointers to this process' graph entries.
size_t getLocalNumVertices() const
Returns the number of vertices on this process.