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 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
223
224 {
225 // Do the import with the Tpetra::CrsMatrix traits object
226 const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
227 RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
228
229 RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
230 *tm, numLocalRows, myNewRows);
231
232 RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
233
234 return xmnew;
235 }
236 }
237};
238
239
241// Tpetra::CrsGraph
242template <typename lno_t,
243 typename gno_t,
244 typename node_t>
245struct XpetraTraits<Tpetra::CrsGraph<lno_t, gno_t, node_t> >
246{
247 typedef typename Xpetra::CrsGraph<lno_t, gno_t, node_t> xgraph_t;
248 typedef typename Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xtgraph_t;
249 typedef typename Tpetra::CrsGraph<lno_t, gno_t, node_t> tgraph_t;
250
251 static inline RCP<xgraph_t> convertToXpetra(const RCP<tgraph_t> &a)
252 {
253 return rcp(new xtgraph_t(a));
254 }
255
256 static RCP<tgraph_t> doMigration(const tgraph_t &from,
257 size_t numLocalRows, const gno_t *myNewRows)
258 {
259 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
260
261 // source map
262 const RCP<const map_t> &smap = from.getRowMap();
263 int oldNumElts = smap->getLocalNumElements();
264 gno_t numGlobalRows = smap->getGlobalNumElements();
265 gno_t base = smap->getMinAllGlobalIndex();
266
267 // target map
268 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
269 const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
270 RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
271
272 // importer
273 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
274
275 // number of entries in my new rows
276 typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
277 vector_t numOld(smap);
278 vector_t numNew(tmap);
279 for (int lid=0; lid < oldNumElts; lid++){
280 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
281 from.getNumEntriesInLocalRow(lid));
282 }
283 numNew.doImport(numOld, importer, Tpetra::INSERT);
284
285 size_t numElts = tmap->getLocalNumElements();
286 ArrayRCP<const gno_t> nnz;
287 if (numElts > 0)
288 nnz = numNew.getData(0); // hangs if vector len == 0
289
290 ArrayRCP<const size_t> nnz_size_t;
291
292 if (numElts && sizeof(gno_t) != sizeof(size_t)){
293 size_t *vals = new size_t [numElts];
294 nnz_size_t = arcp(vals, 0, numElts, true);
295 for (size_t i=0; i < numElts; i++){
296 vals[i] = static_cast<size_t>(nnz[i]);
297 }
298 }
299 else{
300 nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
301 }
302
303 // target graph
304 RCP<tgraph_t> G = rcp(new tgraph_t(tmap, nnz_size_t()));
305
306 G->doImport(from, importer, Tpetra::INSERT);
307 G->fillComplete();
308
309 return G;
310 }
311
312};
313
314
316// Xpetra::CrsGraph
317template <typename lno_t,
318 typename gno_t,
319 typename node_t>
320struct XpetraTraits<Xpetra::CrsGraph<lno_t, gno_t, node_t> >
321{
322 typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
323 typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
324 typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
325
326 static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
327 {
328 return a;
329 }
330
331 static RCP<x_graph_t> doMigration(const x_graph_t &from,
332 size_t numLocalRows, const gno_t *myNewRows)
333 {
334 {
335 // Do the import with the Tpetra::CrsGraph traits object
336 const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
337 RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
338
339 RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
340 *tg, numLocalRows, myNewRows);
341
342 RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
343 return xgnew;
344 }
345 }
346};
347
348
350// Xpetra::CrsGraph specialization
351template < typename node_t>
352struct XpetraTraits<Xpetra::CrsGraph<int, int, node_t> >
353{
354 typedef int lno_t;
355 typedef int gno_t;
356 typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
357 typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
358 typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
359
360 static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
361 {
362 return a;
363 }
364
365 static RCP<x_graph_t> doMigration(const x_graph_t &from,
366 size_t numLocalRows, const gno_t *myNewRows)
367 {
368 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
369
370 {
371 // Do the import with the Tpetra::CrsGraph traits object
372 const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
373 RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
374
375 RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
376 *tg, numLocalRows, myNewRows);
377
378 RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
379
380 return xgnew;
381 }
382 }
383};
384
386// Tpetra::Vector
387template <typename scalar_t,
388 typename lno_t,
389 typename gno_t,
390 typename node_t>
391struct XpetraTraits<Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
392{
393 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
394 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
395 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
396
397 static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
398 {
399 return rcp(new xt_vector_t(a));
400 }
401
402 static RCP<t_vector_t> doMigration(const t_vector_t &from,
403 size_t numLocalElts, const gno_t *myNewElts)
404 {
405 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
406
407 // source map
408 const RCP<const map_t> &smap = from.getMap();
409 gno_t numGlobalElts = smap->getGlobalNumElements();
410 gno_t base = smap->getMinAllGlobalIndex();
411
412 // target map
413 ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
414 const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
415 RCP<const map_t> tmap = rcp(new map_t(numGlobalElts, eltList, base, comm));
416
417 // importer
418 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
419
420 // target vector
421 RCP<t_vector_t> V =
422 Tpetra::createVector<scalar_t,lno_t,gno_t,node_t>(tmap);
423 V->doImport(from, importer, Tpetra::INSERT);
424
425 return V;
426 }
427};
428
429
431// Xpetra::Vector
432template <typename scalar_t,
433 typename lno_t,
434 typename gno_t,
435 typename node_t>
436struct XpetraTraits<Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
437{
438 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
439 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
440 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
441
442 static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
443 {
444 return a;
445 }
446
447 static RCP<x_vector_t> doMigration(const x_vector_t &from,
448 size_t numLocalRows, const gno_t *myNewRows)
449 {
450 {
451 // Do the import with the Tpetra::Vector traits object
452 const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
453 RCP<const t_vector_t> tv = xtv->getTpetra_Vector();
454
455 RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
456 *tv, numLocalRows, myNewRows);
457
458 RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
459
460 return xvnew;
461 }
462 }
463};
464
466// Xpetra::Vector specialization
467template <typename node_t>
468struct XpetraTraits<Xpetra::Vector<double, int, int, node_t> >
469{
470 typedef double scalar_t;
471 typedef int lno_t;
472 typedef int gno_t;
473 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
474 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
475 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
476
477 static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
478 {
479 return a;
480 }
481
482 static RCP<x_vector_t> doMigration(const x_vector_t &from,
483 size_t numLocalRows, const gno_t *myNewRows)
484 {
485 Xpetra::UnderlyingLib lib = from.getMap()->lib();
486
487 {
488 // Do the import with the Tpetra::Vector traits object
489 const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
490 RCP<t_vector_t> tv = xtv->getTpetra_Vector();
491
492 RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
493 *tv, numLocalRows, myNewRows);
494
495 RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
496
497 return xvnew;
498 }
499 }
500};
501
503// Tpetra::MultiVector
504template <typename scalar_t,
505 typename lno_t,
506 typename gno_t,
507 typename node_t>
508struct XpetraTraits<Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
509{
510 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
511 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
512 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
513
514 static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
515 {
516 return rcp(new xt_vector_t(a));
517 }
518
519 static RCP<t_vector_t> doMigration(const t_vector_t &from,
520 size_t numLocalElts, const gno_t *myNewElts)
521 {
522 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
523
524 // source map
525 const RCP<const map_t> &smap = from.getMap();
526 gno_t numGlobalElts = smap->getGlobalNumElements();
527 gno_t base = smap->getMinAllGlobalIndex();
528
529 // target map
530 ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
531 const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
532 RCP<const map_t> tmap = rcp(new map_t(numGlobalElts, eltList, base, comm));
533
534 // importer
535 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
536
537 // target vector
538 RCP<t_vector_t> MV = rcp(
539 new t_vector_t(tmap, from.getNumVectors(), true));
540 MV->doImport(from, importer, Tpetra::INSERT);
541
542 return MV;
543 }
544};
545
546
548// Xpetra::MultiVector
549template <typename scalar_t,
550 typename lno_t,
551 typename gno_t,
552 typename node_t>
553struct XpetraTraits<Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
554{
555 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
556 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
557 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
558
559 static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
560 {
561 return a;
562 }
563
564 static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
565 size_t numLocalRows, const gno_t *myNewRows)
566 {
567 {
568 // Do the import with the Tpetra::MultiVector traits object
569 const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
570 RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
571
572 RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
573 *tv, numLocalRows, myNewRows);
574
575 RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
576
577 return xvnew;
578 }
579 }
580};
581
583// Xpetra::MultiVector specialization
584template <typename node_t>
585struct XpetraTraits<Xpetra::MultiVector<double, int, int, node_t> >
586{
587 typedef double scalar_t;
588 typedef int lno_t;
589 typedef int gno_t;
590 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
591 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
592 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
593
594 static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
595 {
596 return a;
597 }
598
599 static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
600 size_t numLocalRows, const gno_t *myNewRows)
601 {
602 Xpetra::UnderlyingLib lib = from.getMap()->lib();
603
604 {
605 // Do the import with the Tpetra::MultiVector traits object
606 const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
607 RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
608
609 RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
610 *tv, numLocalRows, myNewRows);
611
612 RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
613
614 return xvnew;
615 }
616 }
617};
618
619#endif // DOXYGEN_SHOULD_SKIP_THIS
620
621} //namespace Zoltan2
622
623#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.