50 #ifndef _ZOLTAN2_XPETRATRAITS_HPP_ 51 #define _ZOLTAN2_XPETRATRAITS_HPP_ 56 #include <Xpetra_EpetraCrsMatrix.hpp> 57 #include <Xpetra_TpetraCrsMatrix.hpp> 58 #include <Xpetra_TpetraRowMatrix.hpp> 59 #include <Xpetra_EpetraVector.hpp> 60 #include <Xpetra_TpetraVector.hpp> 61 #include <Xpetra_EpetraUtils.hpp> 62 #include <Tpetra_Vector.hpp> 91 template <
typename User>
115 size_t numLocalRows,
const gno_t *myNewRows)
117 return Teuchos::null;
121 #ifndef DOXYGEN_SHOULD_SKIP_THIS 125 template <
typename scalar_t,
129 struct XpetraTraits<Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
131 typedef typename Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t>
133 typedef typename Xpetra::TpetraCrsMatrix<scalar_t,lno_t,gno_t,node_t>
135 typedef typename Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t>
140 return rcp(
new xtmatrix_t(a));
143 static RCP<tmatrix_t>
doMigration(
const tmatrix_t &from,
144 size_t numLocalRows,
const gno_t *myNewRows)
146 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
150 const RCP<const map_t> &smap = from.getRowMap();
151 int oldNumElts = smap->getNodeNumElements();
152 gno_t numGlobalRows = smap->getGlobalNumElements();
155 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
156 const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
157 RCP<const map_t> tmap = rcp(
158 new map_t(numGlobalRows, rowList, base, comm));
159 int newNumElts = numLocalRows;
162 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
165 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
166 vector_t numOld(smap);
167 vector_t numNew(tmap);
168 for (
int lid=0; lid < oldNumElts; lid++){
169 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
170 scalar_t(from.getNumEntriesInLocalRow(lid)));
172 numNew.doImport(numOld, importer, Tpetra::INSERT);
175 ArrayRCP<size_t> nnz(newNumElts);
177 ArrayRCP<scalar_t> ptr = numNew.getDataNonConst(0);
178 for (
int lid=0; lid < newNumElts; lid++){
179 nnz[lid] =
static_cast<size_t>(ptr[lid]);
184 RCP<tmatrix_t> M = rcp(
new tmatrix_t(tmap, nnz, Tpetra::StaticProfile));
185 M->doImport(from, importer, Tpetra::INSERT);
203 static inline RCP<Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> >
206 return rcp(
new Xpetra::EpetraCrsMatrix(a));
210 static RCP<Epetra_CrsMatrix>
doMigration(
const Epetra_CrsMatrix &from,
211 size_t numLocalRows,
const gno_t *myNewRows)
216 const Epetra_Map &smap = from.RowMap();
217 int oldNumElts = smap.NumMyElements();
218 gno_t numGlobalRows = smap.NumGlobalElements();
221 const Epetra_Comm &comm = from.Comm();
222 Epetra_Map tmap(numGlobalRows, numLocalRows, myNewRows, base, comm);
223 int newNumElts = numLocalRows;
226 Epetra_Import importer(tmap, smap);
229 Epetra_Vector numOld(smap);
230 Epetra_Vector numNew(tmap);
232 for (
int lid=0; lid < oldNumElts; lid++){
233 numOld[lid] = from.NumMyEntries(lid);
235 numNew.Import(numOld, importer, Insert);
237 Array<int> nnz(newNumElts);
238 for (
int lid=0; lid < newNumElts; lid++){
239 nnz[lid] =
static_cast<int>(numNew[lid]);
243 RCP<Epetra_CrsMatrix> M = rcp(
new Epetra_CrsMatrix(::Copy, tmap, nnz.getRawPtr(),
true));
244 M->Import(from, importer, Insert);
255 template <
typename scalar_t,
259 struct XpetraTraits<Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
261 typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
262 typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
263 typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
265 static inline RCP<x_matrix_t>
convertToXpetra(
const RCP<x_matrix_t > &a)
270 static RCP<x_matrix_t>
doMigration(
const x_matrix_t &from,
271 size_t numLocalRows,
const gno_t *myNewRows)
273 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
275 if (lib == Xpetra::UseEpetra){
276 throw std::logic_error(
"compiler should have used specialization");
279 const xt_matrix_t *xtm =
dynamic_cast<const xt_matrix_t *
>(&from);
280 RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
283 *tm, numLocalRows, myNewRows);
295 template <
typename node_t>
296 struct XpetraTraits<Xpetra::CrsMatrix<double, int, int, node_t> >
298 typedef double scalar_t;
301 typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
302 typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
303 typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
304 typedef Xpetra::EpetraCrsMatrix xe_matrix_t;
305 typedef Epetra_CrsMatrix e_matrix_t;
307 static inline RCP<x_matrix_t>
convertToXpetra(
const RCP<x_matrix_t > &a)
312 static RCP<x_matrix_t>
doMigration(
const x_matrix_t &from,
313 size_t numLocalRows,
const gno_t *myNewRows)
315 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
317 if (lib == Xpetra::UseEpetra){
319 const xe_matrix_t *xem =
dynamic_cast<const xe_matrix_t *
>(&from);
320 RCP<const e_matrix_t> em = xem->getEpetra_CrsMatrix();
323 *em, numLocalRows, myNewRows);
330 const xt_matrix_t *xtm =
dynamic_cast<const xt_matrix_t *
>(&from);
331 RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
334 *tm, numLocalRows, myNewRows);
346 template <
typename lno_t,
349 struct XpetraTraits<Tpetra::CrsGraph<lno_t, gno_t, node_t> >
351 typedef typename Xpetra::CrsGraph<lno_t, gno_t, node_t>
xgraph_t;
352 typedef typename Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xtgraph_t;
353 typedef typename Tpetra::CrsGraph<lno_t, gno_t, node_t>
tgraph_t;
357 return rcp(
new xtgraph_t(a));
360 static RCP<tgraph_t>
doMigration(
const tgraph_t &from,
361 size_t numLocalRows,
const gno_t *myNewRows)
363 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
367 const RCP<const map_t> &smap = from.getRowMap();
368 int oldNumElts = smap->getNodeNumElements();
369 gno_t numGlobalRows = smap->getGlobalNumElements();
372 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
373 const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
374 RCP<const map_t> tmap = rcp(
375 new map_t(numGlobalRows, rowList, base, comm));
378 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
381 typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
382 vector_t numOld(smap);
383 vector_t numNew(tmap);
384 for (
int lid=0; lid < oldNumElts; lid++){
385 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
386 from.getNumEntriesInLocalRow(lid));
388 numNew.doImport(numOld, importer, Tpetra::INSERT);
390 size_t numElts = tmap->getNodeNumElements();
391 ArrayRCP<const gno_t> nnz;
393 nnz = numNew.getData(0);
395 ArrayRCP<const size_t> nnz_size_t;
397 if (numElts &&
sizeof(gno_t) !=
sizeof(
size_t)){
398 size_t *
vals =
new size_t [numElts];
399 nnz_size_t = arcp(vals, 0, numElts,
true);
400 for (
size_t i=0; i < numElts; i++){
401 vals[i] =
static_cast<size_t>(nnz[i]);
405 nnz_size_t = arcp_reinterpret_cast<
const size_t>(nnz);
409 RCP<tgraph_t> G = rcp(
new tgraph_t(tmap, nnz_size_t, Tpetra::StaticProfile));
411 G->doImport(from, importer, Tpetra::INSERT);
427 static inline RCP<Xpetra::CrsGraph<lno_t,gno_t,node_t> >
430 return rcp(
new Xpetra::EpetraCrsGraph(a));
433 static RCP<Epetra_CrsGraph>
doMigration(
const Epetra_CrsGraph &from,
434 size_t numLocalRows,
const gno_t *myNewRows)
439 const Epetra_BlockMap &smap = from.RowMap();
440 gno_t numGlobalRows = smap.NumGlobalElements();
441 lno_t oldNumElts = smap.NumMyElements();
444 const Epetra_Comm &comm = from.Comm();
445 Epetra_BlockMap tmap(numGlobalRows, numLocalRows,
446 myNewRows, 1, base, comm);
447 lno_t newNumElts = tmap.NumMyElements();
450 Epetra_Import importer(tmap, smap);
453 Epetra_Vector numOld(smap);
454 Epetra_Vector numNew(tmap);
456 for (
int lid=0; lid < oldNumElts; lid++){
457 numOld[lid] = from.NumMyIndices(lid);
459 numNew.Import(numOld, importer, Insert);
461 Array<int> nnz(newNumElts);
462 for (
int lid=0; lid < newNumElts; lid++){
463 nnz[lid] =
static_cast<int>(numNew[lid]);
467 RCP<Epetra_CrsGraph> G = rcp(
new Epetra_CrsGraph(::Copy, tmap, nnz.getRawPtr(),
true));
468 G->Import(from, importer, Insert);
480 template <
typename lno_t,
483 struct XpetraTraits<Xpetra::CrsGraph<lno_t, gno_t, node_t> >
485 typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
486 typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
487 typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
494 static RCP<x_graph_t>
doMigration(
const x_graph_t &from,
495 size_t numLocalRows,
const gno_t *myNewRows)
497 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
499 if (lib == Xpetra::UseEpetra){
500 throw std::logic_error(
"compiler should have used specialization");
503 const xt_graph_t *xtg =
dynamic_cast<const xt_graph_t *
>(&from);
504 RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
507 *tg, numLocalRows, myNewRows);
519 template <
typename scalar_t,
523 struct XpetraTraits<Xpetra::RowMatrix<scalar_t, lno_t, gno_t, node_t> >
525 typedef Xpetra::RowMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
526 typedef Xpetra::TpetraRowMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
527 typedef Tpetra::RowMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
529 static inline RCP<x_matrix_t>
convertToXpetra(
const RCP<x_matrix_t > &a)
534 static RCP<x_matrix_t>
doMigration(
const x_matrix_t &from,
535 size_t numLocalRows,
const gno_t *myNewRows)
537 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
539 if (lib == Xpetra::UseEpetra){
540 throw std::logic_error(
"compiler should have used specialization");
543 const xt_matrix_t *xtm =
dynamic_cast<const xt_matrix_t *
>(&from);
544 RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
547 *tm, numLocalRows, myNewRows);
558 template <
typename node_t>
559 struct XpetraTraits<Xpetra::CrsGraph<int, int, node_t> >
563 typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
564 typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
565 typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
566 typedef Xpetra::EpetraCrsGraph xe_graph_t;
567 typedef Epetra_CrsGraph e_graph_t;
574 static RCP<x_graph_t>
doMigration(
const x_graph_t &from,
575 size_t numLocalRows,
const gno_t *myNewRows)
577 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
579 if (lib == Xpetra::UseEpetra){
581 const xe_graph_t *xeg =
dynamic_cast<const xe_graph_t *
>(&from);
582 RCP<const e_graph_t> eg = xeg->getEpetra_CrsGraph();
585 *eg, numLocalRows, myNewRows);
592 const xt_graph_t *xtg =
dynamic_cast<const xt_graph_t *
>(&from);
593 RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
596 *tg, numLocalRows, myNewRows);
607 template <
typename scalar_t,
613 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
614 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
615 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
617 static inline RCP<x_vector_t>
convertToXpetra(
const RCP<t_vector_t> &a)
619 return rcp(
new xt_vector_t(a));
622 static RCP<t_vector_t>
doMigration(
const t_vector_t &from,
623 size_t numLocalElts,
const gno_t *myNewElts)
626 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
629 const RCP<const map_t> &smap = from.getMap();
630 gno_t numGlobalElts = smap->getGlobalNumElements();
633 ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
634 const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
635 RCP<const map_t> tmap = rcp(
636 new map_t(numGlobalElts, eltList, base, comm));
639 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
643 Tpetra::createVector<scalar_t,lno_t,gno_t,node_t>(tmap);
644 V->doImport(from, importer, Tpetra::INSERT);
660 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
662 static inline RCP<x_vector_t>
convertToXpetra(
const RCP<Epetra_Vector> &a)
664 RCP<Xpetra::EpetraVector> xev = rcp(
new Xpetra::EpetraVector(a));
665 return rcp_implicit_cast<x_vector_t>(xev);
668 static RCP<Epetra_Vector>
doMigration(
const Epetra_Vector &from,
669 size_t numLocalElts,
const gno_t *myNewElts)
673 const Epetra_BlockMap &smap = from.Map();
674 gno_t numGlobalElts = smap.NumGlobalElements();
677 const Epetra_Comm &comm = from.Comm();
678 const Epetra_BlockMap tmap(numGlobalElts, numLocalElts, myNewElts,
682 Epetra_Import importer(tmap, smap);
685 RCP<Epetra_Vector> V = rcp(
new Epetra_Vector(tmap,
true));
686 Epetra_CombineMode c = Insert;
687 V->Import(from, importer, c);
695 template <
typename scalar_t,
701 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
702 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
703 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
705 static inline RCP<x_vector_t>
convertToXpetra(
const RCP<x_vector_t> &a)
710 static RCP<x_vector_t>
doMigration(
const x_vector_t &from,
711 size_t numLocalRows,
const gno_t *myNewRows)
713 Xpetra::UnderlyingLib lib = from.getMap()->lib();
715 if (lib == Xpetra::UseEpetra){
716 throw std::logic_error(
"compiler should have used specialization");
719 const xt_vector_t *xtv =
dynamic_cast<const xt_vector_t *
>(&from);
720 RCP<const t_vector_t> tv = xtv->getTpetra_Vector();
723 *tv, numLocalRows, myNewRows);
734 template <
typename node_t>
737 typedef double scalar_t;
740 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
741 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
742 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
743 typedef Xpetra::EpetraVector xe_vector_t;
744 typedef Epetra_Vector e_vector_t;
746 static inline RCP<x_vector_t>
convertToXpetra(
const RCP<x_vector_t> &a)
751 static RCP<x_vector_t>
doMigration(
const x_vector_t &from,
752 size_t numLocalRows,
const gno_t *myNewRows)
754 Xpetra::UnderlyingLib lib = from.getMap()->lib();
756 if (lib == Xpetra::UseEpetra){
758 const xe_vector_t *xev =
dynamic_cast<const xe_vector_t *
>(&from);
759 RCP<const e_vector_t> ev = rcp(xev->getEpetra_Vector());
762 *ev, numLocalRows, myNewRows);
769 const xt_vector_t *xtv =
dynamic_cast<const xt_vector_t *
>(&from);
770 RCP<t_vector_t> tv = xtv->getTpetra_Vector();
773 *tv, numLocalRows, myNewRows);
784 template <
typename scalar_t,
788 struct XpetraTraits<Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
790 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
791 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
792 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
794 static inline RCP<x_vector_t>
convertToXpetra(
const RCP<t_vector_t> &a)
796 return rcp(
new xt_vector_t(a));
799 static RCP<t_vector_t>
doMigration(
const t_vector_t &from,
800 size_t numLocalElts,
const gno_t *myNewElts)
802 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
806 const RCP<const map_t> &smap = from.getMap();
807 gno_t numGlobalElts = smap->getGlobalNumElements();
810 ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
811 const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
812 RCP<const map_t> tmap = rcp(
813 new map_t(numGlobalElts, eltList, base, comm));
816 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
819 RCP<t_vector_t> MV = rcp(
820 new t_vector_t(tmap, from.getNumVectors(),
true));
821 MV->doImport(from, importer, Tpetra::INSERT);
836 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
839 const RCP<Epetra_MultiVector> &a)
841 RCP<Xpetra::EpetraMultiVector> xemv = rcp(
new Xpetra::EpetraMultiVector(a));
842 return rcp_implicit_cast<x_mvector_t>(xemv);
845 static RCP<Epetra_MultiVector>
doMigration(
const Epetra_MultiVector &from,
846 size_t numLocalElts,
const gno_t *myNewElts)
850 const Epetra_BlockMap &smap = from.Map();
851 gno_t numGlobalElts = smap.NumGlobalElements();
854 const Epetra_Comm &comm = from.Comm();
855 const Epetra_BlockMap tmap(numGlobalElts, numLocalElts, myNewElts,
859 Epetra_Import importer(tmap, smap);
862 RCP<Epetra_MultiVector> MV = rcp(
863 new Epetra_MultiVector(tmap, from.NumVectors(),
true));
864 Epetra_CombineMode c = Insert;
865 MV->Import(from, importer, c);
873 template <
typename scalar_t,
877 struct XpetraTraits<Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
879 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
880 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
881 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
883 static inline RCP<x_mvector_t>
convertToXpetra(
const RCP<x_mvector_t> &a)
888 static RCP<x_mvector_t>
doMigration(
const x_mvector_t &from,
889 size_t numLocalRows,
const gno_t *myNewRows)
891 Xpetra::UnderlyingLib lib = from.getMap()->lib();
893 if (lib == Xpetra::UseEpetra){
894 throw std::logic_error(
"compiler should have used specialization");
897 const xt_mvector_t *xtv =
dynamic_cast<const xt_mvector_t *
>(&from);
898 RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
901 *tv, numLocalRows, myNewRows);
912 template <
typename node_t>
913 struct XpetraTraits<Xpetra::MultiVector<double, int, int, node_t> >
915 typedef double scalar_t;
918 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
919 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
920 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
921 typedef Xpetra::EpetraMultiVector xe_mvector_t;
922 typedef Epetra_MultiVector e_mvector_t;
924 static inline RCP<x_mvector_t>
convertToXpetra(
const RCP<x_mvector_t> &a)
929 static RCP<x_mvector_t>
doMigration(
const x_mvector_t &from,
930 size_t numLocalRows,
const gno_t *myNewRows)
932 Xpetra::UnderlyingLib lib = from.getMap()->lib();
934 if (lib == Xpetra::UseEpetra){
936 const xe_mvector_t *xev =
dynamic_cast<const xe_mvector_t *
>(&from);
937 RCP<e_mvector_t> ev = xev->getEpetra_MultiVector();
940 *ev, numLocalRows, myNewRows);
948 const xt_mvector_t *xtv =
dynamic_cast<const xt_mvector_t *
>(&from);
949 RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
952 *tv, numLocalRows, myNewRows);
961 #endif // DOXYGEN_SHOULD_SKIP_THIS 965 #endif // _ZOLTAN2_XPETRATRAITS_HPP_ Defines the traits required for Tpetra, Eptra and Xpetra objects.
default_gno_t gno_t
The objects global ordinal data type.
static RCP< User > doMigration(const User &from, size_t numLocalRows, const gno_t *myNewRows)
Migrate the object Given a user object and a new row distribution, create and return a new user objec...
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
default_lno_t lno_t
The objects local ordinal data type.
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
Gathering definitions used in software development.