57 #include <Teuchos_DefaultComm.hpp>
58 #include <Teuchos_RCP.hpp>
59 #include <Teuchos_Array.hpp>
60 #include <Teuchos_ArrayRCP.hpp>
61 #include <Teuchos_Comm.hpp>
62 #include <Teuchos_VerboseObject.hpp>
63 #include <Tpetra_CrsMatrix.hpp>
64 #include <Tpetra_Vector.hpp>
66 #ifdef HAVE_ZOLTAN_EPETRA
67 #include <Xpetra_EpetraUtils.hpp>
72 using Teuchos::ArrayRCP;
73 using Teuchos::ArrayView;
86 if (nglobalrows !=
size_t(maxgid - basegid + 1)){
87 std::cout <<
"Error: Map is invalid for test - fix test" << std::endl;
88 std::cerr <<
"Error: Map is invalid for test - fix test" << std::endl;
89 std::cout <<
"FAIL" << std::endl;
92 RCP<Array<zgno_t> > mygids = rcp(
new Array<zgno_t>);
94 if (firstzgno_t < basegid){
97 firstzgno_t = basegid - n + proc;
99 firstzgno_t = basegid;
101 for (
zgno_t gid=firstzgno_t; gid <= maxgid; gid+=nprocs){
102 (*mygids).append(gid);
105 ArrayRCP<zgno_t> newIdArcp = Teuchos::arcp(mygids);
110 #ifdef HAVE_EPETRA_DATA_TYPES
113 const Epetra_Comm &comm = emap.Comm();
114 int proc = comm.MyPID();
115 int nprocs = comm.NumProc();
116 zgno_t basegid = emap.MinAllGID();
117 zgno_t maxgid = emap.MaxAllGID();
118 size_t nglobalrows = emap.NumGlobalElements();
124 ArrayRCP<zgno_t>
roundRobinMap(
const Tpetra::Map<zlno_t, zgno_t, znode_t> &tmap)
126 const RCP<const Comm<int> > &comm = tmap.getComm();
127 int proc = comm->getRank();
128 int nprocs = comm->getSize();
129 zgno_t basegid = tmap.getMinAllGlobalIndex();
130 zgno_t maxgid = tmap.getMaxAllGlobalIndex();
131 size_t nglobalrows = tmap.getGlobalNumElements();
136 ArrayRCP<zgno_t>
roundRobinMap(
const Xpetra::Map<zlno_t, zgno_t, znode_t> &xmap)
138 const RCP<const Comm<int> > &comm = xmap.getComm();
139 int proc = comm->getRank();
140 int nprocs = comm->getSize();
141 zgno_t basegid = xmap.getMinAllGlobalIndex();
142 zgno_t maxgid = xmap.getMaxAllGlobalIndex();
143 size_t nglobalrows = xmap.getGlobalNumElements();
148 int main(
int narg,
char *arg[])
150 Tpetra::ScopeGuard tscope(&narg, &arg);
151 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
153 int rank = comm->getRank();
156 Teuchos::RCP<Teuchos::FancyOStream> outStream =
157 Teuchos::VerboseObjectBase::getDefaultOStream();
158 Teuchos::EVerbosityLevel v=Teuchos::VERB_EXTREME;
160 typedef Tpetra::CrsMatrix<zscalar_t,zlno_t,zgno_t,znode_t>
tmatrix_t;
161 typedef Tpetra::CrsGraph<zlno_t,zgno_t,znode_t>
tgraph_t;
162 typedef Tpetra::Vector<zscalar_t,zlno_t,zgno_t,znode_t>
tvector_t;
163 typedef Tpetra::MultiVector<zscalar_t,zlno_t,zgno_t,znode_t> tmvector_t;
164 typedef Xpetra::CrsMatrix<zscalar_t,zlno_t,zgno_t,znode_t>
xmatrix_t;
165 typedef Xpetra::CrsGraph<zlno_t,zgno_t,znode_t>
xgraph_t;
166 typedef Xpetra::Vector<zscalar_t,zlno_t,zgno_t,znode_t>
xvector_t;
167 typedef Xpetra::MultiVector<zscalar_t,zlno_t,zgno_t,znode_t> xmvector_t;
171 RCP<UserInputForTests> uinput;
177 catch(std::exception &e){
179 std::cout << e.what() << std::endl;
196 M = uinput->getUITpetraCrsMatrix();
198 catch(std::exception &e){
200 std::cout << e.what() << std::endl;
205 std::cout <<
"Original Tpetra matrix " << M->getGlobalNumRows()
206 <<
" x " << M->getGlobalNumCols() << std::endl;
208 M->describe(*outStream,v);
210 ArrayRCP<zgno_t> newRowIds =
roundRobinMap(*(M->getRowMap()));
212 zgno_t localNumRows = newRowIds.size();
214 RCP<const tmatrix_t> newM;
217 localNumRows, newRowIds.getRawPtr());
219 catch(std::exception &e){
221 std::cout << e.what() << std::endl;
224 " Zoltan2::XpetraTraits<tmatrix_t>::doMigration ", 1);
227 std::cout <<
"Migrated Tpetra matrix" << std::endl;
229 newM->describe(*outStream,v);
237 G = uinput->getUITpetraCrsGraph();
239 catch(std::exception &e){
241 std::cout << e.what() << std::endl;
246 std::cout <<
"Original Tpetra graph" << std::endl;
248 G->describe(*outStream,v);
250 ArrayRCP<zgno_t> newRowIds =
roundRobinMap(*(G->getRowMap()));
252 zgno_t localNumRows = newRowIds.size();
254 RCP<const tgraph_t> newG;
257 localNumRows, newRowIds.getRawPtr());
259 catch(std::exception &e){
261 std::cout << e.what() << std::endl;
264 " Zoltan2::XpetraTraits<tgraph_t>::doMigration ", 1);
267 std::cout <<
"Migrated Tpetra graph" << std::endl;
269 newG->describe(*outStream,v);
277 V = rcp(
new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 1));
280 catch(std::exception &e){
282 std::cout << e.what() << std::endl;
287 std::cout <<
"Original Tpetra vector" << std::endl;
289 V->describe(*outStream,v);
293 zgno_t localNumRows = newRowIds.size();
295 RCP<const tvector_t> newV;
298 localNumRows, newRowIds.getRawPtr());
300 catch(std::exception &e){
302 std::cout << e.what() << std::endl;
305 " Zoltan2::XpetraTraits<tvector_t>::doMigration ", 1);
308 std::cout <<
"Migrated Tpetra vector" << std::endl;
310 newV->describe(*outStream,v);
318 MV = rcp(
new tmvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 3));
321 catch(std::exception &e){
323 std::cout << e.what() << std::endl;
328 std::cout <<
"Original Tpetra multivector" << std::endl;
330 MV->describe(*outStream,v);
334 zgno_t localNumRows = newRowIds.size();
336 RCP<const tmvector_t> newMV;
339 localNumRows, newRowIds.getRawPtr());
341 catch(std::exception &e){
343 std::cout << e.what() << std::endl;
346 " Zoltan2::XpetraTraits<tmvector_t>::doMigration ", 1);
349 std::cout <<
"Migrated Tpetra multivector" << std::endl;
351 newMV->describe(*outStream,v);
366 M = uinput->getUIXpetraCrsMatrix();
368 catch(std::exception &e){
370 std::cout << e.what() << std::endl;
375 std::cout <<
"Original Xpetra matrix" << std::endl;
377 M->describe(*outStream,v);
379 ArrayRCP<zgno_t> newRowIds =
roundRobinMap(*(M->getRowMap()));
381 zgno_t localNumRows = newRowIds.size();
383 RCP<const xmatrix_t> newM;
386 localNumRows, newRowIds.getRawPtr());
388 catch(std::exception &e){
390 std::cout << e.what() << std::endl;
393 " Zoltan2::XpetraTraits<xmatrix_t>::doMigration ", 1);
396 std::cout <<
"Migrated Xpetra matrix" << std::endl;
398 newM->describe(*outStream,v);
406 G = uinput->getUIXpetraCrsGraph();
408 catch(std::exception &e){
410 std::cout << e.what() << std::endl;
415 std::cout <<
"Original Xpetra graph" << std::endl;
417 G->describe(*outStream,v);
419 ArrayRCP<zgno_t> newRowIds =
roundRobinMap(*(G->getRowMap()));
421 zgno_t localNumRows = newRowIds.size();
423 RCP<const xgraph_t> newG;
426 localNumRows, newRowIds.getRawPtr());
428 catch(std::exception &e){
430 std::cout << e.what() << std::endl;
433 " Zoltan2::XpetraTraits<xgraph_t>::doMigration ", 1);
436 std::cout <<
"Migrated Xpetra graph" << std::endl;
438 newG->describe(*outStream,v);
447 rcp(
new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 1));
451 catch(std::exception &e){
453 std::cout << e.what() << std::endl;
458 std::cout <<
"Original Xpetra vector" << std::endl;
460 V->describe(*outStream,v);
464 zgno_t localNumRows = newRowIds.size();
466 RCP<const xvector_t> newV;
469 localNumRows, newRowIds.getRawPtr());
471 catch(std::exception &e){
473 std::cout << e.what() << std::endl;
476 " Zoltan2::XpetraTraits<xvector_t>::doMigration ", 1);
479 std::cout <<
"Migrated Xpetra vector" << std::endl;
481 newV->describe(*outStream,v);
489 RCP<tmvector_t> tMV =
490 rcp(
new tmvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 3));
494 catch(std::exception &e){
496 std::cout << e.what() << std::endl;
501 std::cout <<
"Original Xpetra multivector" << std::endl;
503 MV->describe(*outStream,v);
507 zgno_t localNumRows = newRowIds.size();
509 RCP<const xmvector_t> newMV;
512 localNumRows, newRowIds.getRawPtr());
514 catch(std::exception &e){
516 std::cout << e.what() << std::endl;
519 " Zoltan2::XpetraTraits<xmvector_t>::doMigration ", 1);
522 std::cout <<
"Migrated Xpetra multivector" << std::endl;
524 newMV->describe(*outStream,v);
527 #ifdef HAVE_EPETRA_DATA_TYPES
535 typedef Epetra_CrsMatrix ematrix_t;
536 typedef Epetra_CrsGraph egraph_t;
537 typedef Epetra_Vector evector_t;
538 typedef Epetra_MultiVector emvector_t;
539 typedef Epetra_BlockMap emap_t;
543 RCP<UserInputForTests> euinput;
549 catch(std::exception &e){
551 std::cout << e.what() << std::endl;
560 M = euinput->getUIEpetraCrsMatrix();
562 catch(std::exception &e){
564 std::cout << e.what() << std::endl;
569 std::cout <<
"Original Epetra matrix" << std::endl;
573 RCP<const emap_t> emap = Teuchos::rcpFromRef(M->RowMap());
576 zgno_t localNumRows = newRowIds.size();
578 RCP<const ematrix_t> newM;
581 localNumRows, newRowIds.getRawPtr());
583 catch(std::exception &e){
585 std::cout << e.what() << std::endl;
588 " Zoltan2::XpetraTraits<ematrix_t>::doMigration ", 1);
591 std::cout <<
"Migrated Epetra matrix" << std::endl;
593 newM->Print(std::cout);
601 G = euinput->getUIEpetraCrsGraph();
603 catch(std::exception &e){
605 std::cout << e.what() << std::endl;
610 std::cout <<
"Original Epetra graph" << std::endl;
614 RCP<const emap_t> emap = Teuchos::rcpFromRef(G->RowMap());
617 zgno_t localNumRows = newRowIds.size();
619 RCP<const egraph_t> newG;
622 localNumRows, newRowIds.getRawPtr());
624 catch(std::exception &e){
626 std::cout << e.what() << std::endl;
629 " Zoltan2::XpetraTraits<egraph_t>::doMigration ", 1);
632 std::cout <<
"Migrated Epetra graph" << std::endl;
634 newG->Print(std::cout);
642 V = rcp(
new Epetra_Vector(euinput->getUIEpetraCrsGraph()->RowMap()));
645 catch(std::exception &e){
647 std::cout << e.what() << std::endl;
652 std::cout <<
"Original Epetra vector" << std::endl;
656 RCP<const emap_t> emap = Teuchos::rcpFromRef(V->Map());
659 zgno_t localNumRows = newRowIds.size();
661 RCP<const evector_t> newV;
664 localNumRows, newRowIds.getRawPtr());
666 catch(std::exception &e){
668 std::cout << e.what() << std::endl;
671 " Zoltan2::XpetraTraits<evector_t>::doMigration ", 1);
674 std::cout <<
"Migrated Epetra vector" << std::endl;
676 newV->Print(std::cout);
685 rcp(
new Epetra_MultiVector(euinput->getUIEpetraCrsGraph()->RowMap(),3));
688 catch(std::exception &e){
690 std::cout << e.what() << std::endl;
695 std::cout <<
"Original Epetra multivector" << std::endl;
697 MV->Print(std::cout);
699 RCP<const emap_t> emap = Teuchos::rcpFromRef(MV->Map());
702 zgno_t localNumRows = newRowIds.size();
704 RCP<const emvector_t> newMV;
707 localNumRows, newRowIds.getRawPtr());
709 catch(std::exception &e){
711 std::cout << e.what() << std::endl;
714 " Zoltan2::XpetraTraits<emvector_t>::doMigration ", 1);
717 std::cout <<
"Migrated Epetra multivector" << std::endl;
719 newMV->Print(std::cout);
721 #endif // have epetra data types (int, int, double)
728 std::cout <<
"PASS" << std::endl;