Zoltan2
Loading...
Searching...
No Matches
Zoltan2_XpetraTraits.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14#ifndef _ZOLTAN2_XPETRATRAITS_HPP_
15#define _ZOLTAN2_XPETRATRAITS_HPP_
16
18#include <Zoltan2_Standards.hpp>
19
20#include <Xpetra_TpetraCrsGraph.hpp>
21#include <Xpetra_TpetraCrsMatrix.hpp>
22#include <Xpetra_TpetraVector.hpp>
23#include <Tpetra_Vector.hpp>
24
25namespace Zoltan2 {
26
28// Extra traits needed only for Xpetra matrices and graphs
29
48template <typename User>
50{
53 static inline RCP<User> convertToXpetra(const RCP<User> &a)
54 {
55 return a;
56 }
57
61
65
71 static RCP<User> doMigration(const User &from,
72 size_t numLocalRows, const gno_t *myNewRows)
73 {
74 return Teuchos::null;
75 }
76};
77
78#ifndef DOXYGEN_SHOULD_SKIP_THIS
79
81// Tpetra::CrsMatrix
82template <typename scalar_t,
83 typename lno_t,
84 typename gno_t,
85 typename node_t>
86struct XpetraTraits<Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
87{
88 typedef typename Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> xmatrix_t;
89 typedef typename Xpetra::TpetraCrsMatrix<scalar_t,lno_t,gno_t,node_t> xtmatrix_t;
90 typedef typename Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> tmatrix_t;
91
92 static inline RCP<xmatrix_t> convertToXpetra(const RCP<tmatrix_t> &a)
93 {
94 return rcp(new xtmatrix_t(a));
95 }
96
97 static RCP<tmatrix_t> doMigration(const tmatrix_t &from,
98 size_t numLocalRows, const gno_t *myNewRows)
99 {
100 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
101
102 // source map
103 const RCP<const map_t> &smap = from.getRowMap();
104 gno_t numGlobalRows = smap->getGlobalNumElements();
105 gno_t base = smap->getMinAllGlobalIndex();
106
107 // target map
108 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
109 const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
110 RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
111
112 // importer
113 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
114
115 // target matrix
116 // Chris Siefert proposed using the following to make migration
117 // more efficient.
118 // By default, the Domain and Range maps are the same as in "from".
119 // As in the original code, we instead set them both to tmap.
120 // The assumption is a square matrix.
121 // TODO: what about rectangular matrices?
122 // TODO: Should choice of domain/range maps be an option to this function?
123
124 // KDD 3/7/16: disabling Chris' new code to avoid dashboard failures;
125 // KDD 3/7/16: can re-enable when issue #114 is fixed.
126 // KDD 3/7/16: when re-enable CSIEFERT code, can comment out
127 // KDD 3/7/16: "Original way" code.
128 // KDD 1/27/17: Re-enabling Chris' code, as this issue is resolved.
129 RCP<tmatrix_t> M;
130 from.importAndFillComplete(M, importer, tmap, tmap);
131
134 //int oldNumElts = smap->getLocalNumElements();
135 //int newNumElts = numLocalRows;
136
138 //typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
139 //vector_t numOld(smap); // TODO These vectors should have scalar=size_t,
140 //vector_t numNew(tmap); // but ETI does not yet support that.
141 //for (int lid=0; lid < oldNumElts; lid++){
142 //numOld.replaceGlobalValue(smap->getGlobalElement(lid),
143 //scalar_t(from.getNumEntriesInLocalRow(lid)));
144 //}
145 //numNew.doImport(numOld, importer, Tpetra::INSERT);
146
147 // TODO Could skip this copy if could declare vector with scalar=size_t.
148 //ArrayRCP<size_t> nnz(newNumElts);
149 //if (newNumElts > 0){
150 //ArrayRCP<scalar_t> ptr = numNew.getDataNonConst(0);
151 //for (int lid=0; lid < newNumElts; lid++){
152 //nnz[lid] = static_cast<size_t>(ptr[lid]);
153 //}
154 //}
155
156 //RCP<tmatrix_t> M = rcp(new tmatrix_t(tmap, nnz, Tpetra::StaticProfile));
157 //M->doImport(from, importer, Tpetra::INSERT);
158 //M->fillComplete();
159
160 // End of original way we did it.
161 return M;
162 }
163};
164
165
167// Xpetra::CrsMatrix
168template <typename scalar_t,
169 typename lno_t,
170 typename gno_t,
171 typename node_t>
172struct XpetraTraits<Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
173{
174 typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
175 typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
176 typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
177
178 static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
179 {
180 return a;
181 }
182
183 static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
184 size_t numLocalRows, const gno_t *myNewRows)
185 {
186 {
187 // Do the import with the Tpetra::CrsMatrix traits object
188 const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
189 RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
190
191 RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
192 *tm, numLocalRows, myNewRows);
193
194 RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
195
196 return xmnew;
197 }
198 }
199};
200
202// Xpetra::CrsMatrix specialization
203
204template <typename node_t>
205struct XpetraTraits<Xpetra::CrsMatrix<double, int, int, node_t> >
206{
207 typedef double scalar_t;
208 typedef int lno_t;
209 typedef int gno_t;
210 typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
211 typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
212 typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
213
214 static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
215 {
216 return a;
217 }
218
219 static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
220 size_t numLocalRows, const gno_t *myNewRows)
221 {
222 {
223 // Do the import with the Tpetra::CrsMatrix traits object
224 const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
225 RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
226
227 RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
228 *tm, numLocalRows, myNewRows);
229
230 RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
231
232 return xmnew;
233 }
234 }
235};
236
237
239// Tpetra::CrsGraph
240template <typename lno_t,
241 typename gno_t,
242 typename node_t>
243struct XpetraTraits<Tpetra::CrsGraph<lno_t, gno_t, node_t> >
244{
245 typedef typename Xpetra::CrsGraph<lno_t, gno_t, node_t> xgraph_t;
246 typedef typename Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xtgraph_t;
247 typedef typename Tpetra::CrsGraph<lno_t, gno_t, node_t> tgraph_t;
248
249 static inline RCP<xgraph_t> convertToXpetra(const RCP<tgraph_t> &a)
250 {
251 return rcp(new xtgraph_t(a));
252 }
253
254 static RCP<tgraph_t> doMigration(const tgraph_t &from,
255 size_t numLocalRows, const gno_t *myNewRows)
256 {
257 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
258
259 // source map
260 const RCP<const map_t> &smap = from.getRowMap();
261 int oldNumElts = smap->getLocalNumElements();
262 gno_t numGlobalRows = smap->getGlobalNumElements();
263 gno_t base = smap->getMinAllGlobalIndex();
264
265 // target map
266 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
267 const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
268 RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
269
270 // importer
271 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
272
273 // number of entries in my new rows
274 typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
275 vector_t numOld(smap);
276 vector_t numNew(tmap);
277 for (int lid=0; lid < oldNumElts; lid++){
278 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
279 from.getNumEntriesInLocalRow(lid));
280 }
281 numNew.doImport(numOld, importer, Tpetra::INSERT);
282
283 size_t numElts = tmap->getLocalNumElements();
284 ArrayRCP<const gno_t> nnz;
285 if (numElts > 0)
286 nnz = numNew.getData(0); // hangs if vector len == 0
287
288 ArrayRCP<const size_t> nnz_size_t;
289
290 if (numElts && sizeof(gno_t) != sizeof(size_t)){
291 size_t *vals = new size_t [numElts];
292 nnz_size_t = arcp(vals, 0, numElts, true);
293 for (size_t i=0; i < numElts; i++){
294 vals[i] = static_cast<size_t>(nnz[i]);
295 }
296 }
297 else{
298 nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
299 }
300
301 // target graph
302 RCP<tgraph_t> G = rcp(new tgraph_t(tmap, nnz_size_t()));
303
304 G->doImport(from, importer, Tpetra::INSERT);
305 G->fillComplete();
306
307 return G;
308 }
309
310};
311
312
314// Xpetra::CrsGraph
315template <typename lno_t,
316 typename gno_t,
317 typename node_t>
318struct XpetraTraits<Xpetra::CrsGraph<lno_t, gno_t, node_t> >
319{
320 typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
321 typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
322 typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
323
324 static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
325 {
326 return a;
327 }
328
329 static RCP<x_graph_t> doMigration(const x_graph_t &from,
330 size_t numLocalRows, const gno_t *myNewRows)
331 {
332 {
333 // Do the import with the Tpetra::CrsGraph traits object
334 const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
335 RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
336
337 RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
338 *tg, numLocalRows, myNewRows);
339
340 RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
341 return xgnew;
342 }
343 }
344};
345
346
348// Xpetra::CrsGraph specialization
349template < typename node_t>
350struct XpetraTraits<Xpetra::CrsGraph<int, int, node_t> >
351{
352 typedef int lno_t;
353 typedef int gno_t;
354 typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
355 typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
356 typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
357
358 static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
359 {
360 return a;
361 }
362
363 static RCP<x_graph_t> doMigration(const x_graph_t &from,
364 size_t numLocalRows, const gno_t *myNewRows)
365 {
366 {
367 // Do the import with the Tpetra::CrsGraph traits object
368 const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
369 RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
370
371 RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
372 *tg, numLocalRows, myNewRows);
373
374 RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
375
376 return xgnew;
377 }
378 }
379};
380
382// Tpetra::Vector
383template <typename scalar_t,
384 typename lno_t,
385 typename gno_t,
386 typename node_t>
387struct XpetraTraits<Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
388{
389 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
390 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
391 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
392
393 static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
394 {
395 return rcp(new xt_vector_t(a));
396 }
397
398 static RCP<t_vector_t> doMigration(const t_vector_t &from,
399 size_t numLocalElts, const gno_t *myNewElts)
400 {
401 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
402
403 // source map
404 const RCP<const map_t> &smap = from.getMap();
405 gno_t numGlobalElts = smap->getGlobalNumElements();
406 gno_t base = smap->getMinAllGlobalIndex();
407
408 // target map
409 ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
410 const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
411 RCP<const map_t> tmap = rcp(new map_t(numGlobalElts, eltList, base, comm));
412
413 // importer
414 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
415
416 // target vector
417 RCP<t_vector_t> V =
418 Tpetra::createVector<scalar_t,lno_t,gno_t,node_t>(tmap);
419 V->doImport(from, importer, Tpetra::INSERT);
420
421 return V;
422 }
423};
424
425
427// Xpetra::Vector
428template <typename scalar_t,
429 typename lno_t,
430 typename gno_t,
431 typename node_t>
432struct XpetraTraits<Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
433{
434 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
435 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
436 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
437
438 static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
439 {
440 return a;
441 }
442
443 static RCP<x_vector_t> doMigration(const x_vector_t &from,
444 size_t numLocalRows, const gno_t *myNewRows)
445 {
446 {
447 // Do the import with the Tpetra::Vector traits object
448 const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
449 RCP<const t_vector_t> tv = xtv->getTpetra_Vector();
450
451 RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
452 *tv, numLocalRows, myNewRows);
453
454 RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
455
456 return xvnew;
457 }
458 }
459};
460
462// Xpetra::Vector specialization
463template <typename node_t>
464struct XpetraTraits<Xpetra::Vector<double, int, int, node_t> >
465{
466 typedef double scalar_t;
467 typedef int lno_t;
468 typedef int gno_t;
469 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
470 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
471 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
472
473 static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
474 {
475 return a;
476 }
477
478 static RCP<x_vector_t> doMigration(const x_vector_t &from,
479 size_t numLocalRows, const gno_t *myNewRows)
480 {
481 {
482 // Do the import with the Tpetra::Vector traits object
483 const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
484 RCP<t_vector_t> tv = xtv->getTpetra_Vector();
485
486 RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
487 *tv, numLocalRows, myNewRows);
488
489 RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
490
491 return xvnew;
492 }
493 }
494};
495
497// Tpetra::MultiVector
498template <typename scalar_t,
499 typename lno_t,
500 typename gno_t,
501 typename node_t>
502struct XpetraTraits<Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
503{
504 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
505 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
506 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
507
508 static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
509 {
510 return rcp(new xt_vector_t(a));
511 }
512
513 static RCP<t_vector_t> doMigration(const t_vector_t &from,
514 size_t numLocalElts, const gno_t *myNewElts)
515 {
516 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
517
518 // source map
519 const RCP<const map_t> &smap = from.getMap();
520 gno_t numGlobalElts = smap->getGlobalNumElements();
521 gno_t base = smap->getMinAllGlobalIndex();
522
523 // target map
524 ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
525 const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
526 RCP<const map_t> tmap = rcp(new map_t(numGlobalElts, eltList, base, comm));
527
528 // importer
529 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
530
531 // target vector
532 RCP<t_vector_t> MV = rcp(
533 new t_vector_t(tmap, from.getNumVectors(), true));
534 MV->doImport(from, importer, Tpetra::INSERT);
535
536 return MV;
537 }
538};
539
540
542// Xpetra::MultiVector
543template <typename scalar_t,
544 typename lno_t,
545 typename gno_t,
546 typename node_t>
547struct XpetraTraits<Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
548{
549 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
550 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
551 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
552
553 static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
554 {
555 return a;
556 }
557
558 static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
559 size_t numLocalRows, const gno_t *myNewRows)
560 {
561 {
562 // Do the import with the Tpetra::MultiVector traits object
563 const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
564 RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
565
566 RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
567 *tv, numLocalRows, myNewRows);
568
569 RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
570
571 return xvnew;
572 }
573 }
574};
575
577// Xpetra::MultiVector specialization
578template <typename node_t>
579struct XpetraTraits<Xpetra::MultiVector<double, int, int, node_t> >
580{
581 typedef double scalar_t;
582 typedef int lno_t;
583 typedef int gno_t;
584 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
585 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
586 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
587
588 static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
589 {
590 return a;
591 }
592
593 static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
594 size_t numLocalRows, const gno_t *myNewRows)
595 {
596 {
597 // Do the import with the Tpetra::MultiVector traits object
598 const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
599 RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
600
601 RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
602 *tv, numLocalRows, myNewRows);
603
604 RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
605
606 return xvnew;
607 }
608 }
609};
610
611#endif // DOXYGEN_SHOULD_SKIP_THIS
612
613} //namespace Zoltan2
614
615#endif // _ZOLTAN2_XPETRATRAITS_HPP_
typename Zoltan2::InputTraits< ztcrsmatrix_t >::node_t node_t
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tgraph_t
Xpetra::CrsGraph< zlno_t, zgno_t, znode_t > xgraph_t
Xpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > xmatrix_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tmatrix_t
Traits for application input objects.
Gathering definitions used in software development.
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
Definition coloring1.cpp:45
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
Tpetra::Map map_t
Created by mbenlioglu on Aug 31, 2020.
::Tpetra::Details::DefaultTypes::local_ordinal_type default_lno_t
::Tpetra::Details::DefaultTypes::global_ordinal_type default_gno_t
Defines the traits required for Tpetra, Eptra and Xpetra objects.
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_gno_t gno_t
The objects global ordinal data type.
default_lno_t lno_t
The objects local ordinal data type.