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
25#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
26#include <Xpetra_EpetraCrsMatrix.hpp>
27#include <Xpetra_EpetraVector.hpp>
28#include <Xpetra_EpetraUtils.hpp>
29#endif
30
31namespace Zoltan2 {
32
34// Extra traits needed only for Xpetra matrices and graphs
35
58template <typename User>
60{
63 static inline RCP<User> convertToXpetra(const RCP<User> &a)
64 {
65 return a;
66 }
67
71
75
81 static RCP<User> doMigration(const User &from,
82 size_t numLocalRows, const gno_t *myNewRows)
83 {
84 return Teuchos::null;
85 }
86};
87
88#ifndef DOXYGEN_SHOULD_SKIP_THIS
89
91// Tpetra::CrsMatrix
92template <typename scalar_t,
93 typename lno_t,
94 typename gno_t,
95 typename node_t>
96struct XpetraTraits<Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
97{
98 typedef typename Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> xmatrix_t;
99 typedef typename Xpetra::TpetraCrsMatrix<scalar_t,lno_t,gno_t,node_t> xtmatrix_t;
100 typedef typename Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> tmatrix_t;
101
102 static inline RCP<xmatrix_t> convertToXpetra(const RCP<tmatrix_t> &a)
103 {
104 return rcp(new xtmatrix_t(a));
105 }
106
107 static RCP<tmatrix_t> doMigration(const tmatrix_t &from,
108 size_t numLocalRows, const gno_t *myNewRows)
109 {
110 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
111
112 // source map
113 const RCP<const map_t> &smap = from.getRowMap();
114 gno_t numGlobalRows = smap->getGlobalNumElements();
115 gno_t base = smap->getMinAllGlobalIndex();
116
117 // target map
118 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
119 const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
120 RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
121
122 // importer
123 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
124
125 // target matrix
126 // Chris Siefert proposed using the following to make migration
127 // more efficient.
128 // By default, the Domain and Range maps are the same as in "from".
129 // As in the original code, we instead set them both to tmap.
130 // The assumption is a square matrix.
131 // TODO: what about rectangular matrices?
132 // TODO: Should choice of domain/range maps be an option to this function?
133
134 // KDD 3/7/16: disabling Chris' new code to avoid dashboard failures;
135 // KDD 3/7/16: can re-enable when issue #114 is fixed.
136 // KDD 3/7/16: when re-enable CSIEFERT code, can comment out
137 // KDD 3/7/16: "Original way" code.
138 // KDD 1/27/17: Re-enabling Chris' code, as this issue is resolved.
139 RCP<tmatrix_t> M;
140 from.importAndFillComplete(M, importer, tmap, tmap);
141
144 //int oldNumElts = smap->getLocalNumElements();
145 //int newNumElts = numLocalRows;
146
148 //typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
149 //vector_t numOld(smap); // TODO These vectors should have scalar=size_t,
150 //vector_t numNew(tmap); // but ETI does not yet support that.
151 //for (int lid=0; lid < oldNumElts; lid++){
152 //numOld.replaceGlobalValue(smap->getGlobalElement(lid),
153 //scalar_t(from.getNumEntriesInLocalRow(lid)));
154 //}
155 //numNew.doImport(numOld, importer, Tpetra::INSERT);
156
157 // TODO Could skip this copy if could declare vector with scalar=size_t.
158 //ArrayRCP<size_t> nnz(newNumElts);
159 //if (newNumElts > 0){
160 //ArrayRCP<scalar_t> ptr = numNew.getDataNonConst(0);
161 //for (int lid=0; lid < newNumElts; lid++){
162 //nnz[lid] = static_cast<size_t>(ptr[lid]);
163 //}
164 //}
165
166 //RCP<tmatrix_t> M = rcp(new tmatrix_t(tmap, nnz, Tpetra::StaticProfile));
167 //M->doImport(from, importer, Tpetra::INSERT);
168 //M->fillComplete();
169
170 // End of original way we did it.
171 return M;
172 }
173};
174
176#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
177// Epetra_CrsMatrix
178template <>
179struct XpetraTraits<Epetra_CrsMatrix>
180{
185
186 static inline RCP<Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> >
187 convertToXpetra(const RCP<Epetra_CrsMatrix> &a)
188 {
189 RCP<Xpetra::EpetraCrsMatrixT<gno_t, node_t> > xa;
190 try {
191 xa = rcp(new Xpetra::EpetraCrsMatrixT<gno_t, node_t>(a));
192 }
193 catch (std::exception &e) {
194 if (std::is_same<node_t, Xpetra::EpetraNode>::value)
195 throw std::runtime_error(std::string("Cannot convert from "
196 "Epetra_CrsMatrix to "
197 "Xpetra::EpetraCrsMatrixT\n")
198 + e.what());
199 else
200 throw std::runtime_error(std::string("Cannot convert from "
201 "Epetra_CrsMatrix to "
202 "Xpetra::EpetraCrsMatrixT\n"
203 "Use node_t that is supported by "
204 "Xpetra with Epetra classes\n")
205 + e.what());
206 }
207 return xa;
208 }
209
210
211 static RCP<Epetra_CrsMatrix> doMigration(const Epetra_CrsMatrix &from,
212 size_t numLocalRows, const gno_t *myNewRows)
213 {
214 // source map
215 const Epetra_Map &smap = from.RowMap();
216 gno_t numGlobalRows = smap.NumGlobalElements();
217 int base = smap.MinAllGID();
218
219 // target map
220 const Epetra_Comm &comm = from.Comm();
221 Epetra_Map tmap(numGlobalRows, numLocalRows, myNewRows, base, comm);
222
223 // importer
224 Epetra_Import importer(tmap, smap);
225
226
227 // target matrix
228 // Chris Siefert proposed using the following to make migration
229 // more efficient.
230 // By default, the Domain and Range maps are the same as in "from".
231 // As in the original code, we instead set them both to tmap.
232 // The assumption is a square matrix.
233 // TODO: what about rectangular matrices?
234 // TODO: Should choice of domain/range maps be an option to this function?
235
236 RCP<Epetra_CrsMatrix> M = rcp(new Epetra_CrsMatrix(from, importer,
237 &tmap, &tmap));
238
239 // Original way we did it:
240 //
241 // int oldNumElts = smap.NumMyElements();
242 //
243 // // number of non zeros in my new rows
244 // Epetra_Vector numOld(smap);
245 // Epetra_Vector numNew(tmap);
246 //
247 // for (int lid=0; lid < oldNumElts; lid++){
248 // numOld[lid] = from.NumMyEntries(lid);
249 // }
250 // numNew.Import(numOld, importer, Insert);
251 //
252 // int newNumElts = numLocalRows;
253 // Array<int> nnz(newNumElts);
254 // for (int lid=0; lid < newNumElts; lid++){
255 // nnz[lid] = static_cast<int>(numNew[lid]);
256 // }
257 //
258 // RCP<Epetra_CrsMatrix> M =
259 // rcp(new Epetra_CrsMatrix(::Copy, tmap, nnz.getRawPtr(), true));
260 // M->Import(from, importer, Insert);
261 // M->FillComplete();
262
263 return M;
264 }
265};
266#endif
267
269// Xpetra::CrsMatrix
270// KDDKDD: Do we need specializations for Xpetra::EpetraCrsMatrix and
271// KDDKDD: Xpetra::TpetraCrsMatrix
272template <typename scalar_t,
273 typename lno_t,
274 typename gno_t,
275 typename node_t>
276struct XpetraTraits<Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
277{
278 typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
279 typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
280 typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
281
282 static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
283 {
284 return a;
285 }
286
287 static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
288 size_t numLocalRows, const gno_t *myNewRows)
289 {
290 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
291
292 if (lib == Xpetra::UseEpetra){
293 throw std::logic_error("compiler should have used specialization");
294 } else{
295 // Do the import with the Tpetra::CrsMatrix traits object
296 const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
297 RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
298
299 RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
300 *tm, numLocalRows, myNewRows);
301
302 RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
303
304 return xmnew;
305 }
306 }
307};
308
310// Xpetra::CrsMatrix specialization
311
312template <typename node_t>
313struct XpetraTraits<Xpetra::CrsMatrix<double, int, int, node_t> >
314{
315 typedef double scalar_t;
316 typedef int lno_t;
317 typedef int gno_t;
318 typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
319 typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
320 typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
321
322 static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
323 {
324 return a;
325 }
326
327 static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
328 size_t numLocalRows, const gno_t *myNewRows)
329 {
330 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
331
332 if (lib == Xpetra::UseEpetra){
333#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
334 typedef Epetra_CrsMatrix e_matrix_t;
335 typedef Xpetra::EpetraCrsMatrixT<gno_t,node_t> xe_matrix_t;
336 // Do the import with the Epetra_CrsMatrix traits object
337 const xe_matrix_t *xem = dynamic_cast<const xe_matrix_t *>(&from);
338 RCP<const e_matrix_t> em = xem->getEpetra_CrsMatrix();
339
340 RCP<e_matrix_t> emnew = XpetraTraits<e_matrix_t>::doMigration(
341 *em, numLocalRows, myNewRows);
342
343 RCP<x_matrix_t> xmnew = XpetraTraits<e_matrix_t>::convertToXpetra(emnew);
344
345 return xmnew;
346#else
347 throw std::runtime_error("Xpetra with Epetra requested, but "
348 "Trilinos is not built with Epetra");
349#endif
350 } else{
351 // Do the import with the Tpetra::CrsMatrix traits object
352 const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
353 RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
354
355 RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
356 *tm, numLocalRows, myNewRows);
357
358 RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
359
360 return xmnew;
361 }
362 }
363};
364
365
367// Tpetra::CrsGraph
368template <typename lno_t,
369 typename gno_t,
370 typename node_t>
371struct XpetraTraits<Tpetra::CrsGraph<lno_t, gno_t, node_t> >
372{
373 typedef typename Xpetra::CrsGraph<lno_t, gno_t, node_t> xgraph_t;
374 typedef typename Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xtgraph_t;
375 typedef typename Tpetra::CrsGraph<lno_t, gno_t, node_t> tgraph_t;
376
377 static inline RCP<xgraph_t> convertToXpetra(const RCP<tgraph_t> &a)
378 {
379 return rcp(new xtgraph_t(a));
380 }
381
382 static RCP<tgraph_t> doMigration(const tgraph_t &from,
383 size_t numLocalRows, const gno_t *myNewRows)
384 {
385 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
386
387 // source map
388 const RCP<const map_t> &smap = from.getRowMap();
389 int oldNumElts = smap->getLocalNumElements();
390 gno_t numGlobalRows = smap->getGlobalNumElements();
391 gno_t base = smap->getMinAllGlobalIndex();
392
393 // target map
394 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
395 const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
396 RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
397
398 // importer
399 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
400
401 // number of entries in my new rows
402 typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
403 vector_t numOld(smap);
404 vector_t numNew(tmap);
405 for (int lid=0; lid < oldNumElts; lid++){
406 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
407 from.getNumEntriesInLocalRow(lid));
408 }
409 numNew.doImport(numOld, importer, Tpetra::INSERT);
410
411 size_t numElts = tmap->getLocalNumElements();
412 ArrayRCP<const gno_t> nnz;
413 if (numElts > 0)
414 nnz = numNew.getData(0); // hangs if vector len == 0
415
416 ArrayRCP<const size_t> nnz_size_t;
417
418 if (numElts && sizeof(gno_t) != sizeof(size_t)){
419 size_t *vals = new size_t [numElts];
420 nnz_size_t = arcp(vals, 0, numElts, true);
421 for (size_t i=0; i < numElts; i++){
422 vals[i] = static_cast<size_t>(nnz[i]);
423 }
424 }
425 else{
426 nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
427 }
428
429 // target graph
430 RCP<tgraph_t> G = rcp(new tgraph_t(tmap, nnz_size_t()));
431
432 G->doImport(from, importer, Tpetra::INSERT);
433 G->fillComplete();
434
435 return G;
436 }
437
438};
439
441#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
442// Epetra_CrsGraph
443template < >
444struct XpetraTraits<Epetra_CrsGraph>
445{
449 static inline RCP<Xpetra::CrsGraph<lno_t,gno_t,node_t> >
450 convertToXpetra(const RCP<Epetra_CrsGraph> &a)
451 {
452 RCP<Xpetra::EpetraCrsGraphT<gno_t, node_t> > xa;
453 try {
454 xa = rcp(new Xpetra::EpetraCrsGraphT<gno_t, node_t>(a));
455 }
456 catch (std::exception &e) {
457 if (std::is_same<node_t, Xpetra::EpetraNode>::value)
458 throw std::runtime_error(std::string("Cannot convert from "
459 "Epetra_CrsGraph to "
460 "Xpetra::EpetraCrsGraphT\n")
461 + e.what());
462 else
463 throw std::runtime_error(std::string("Cannot convert from "
464 "Epetra_CrsGraph to "
465 "Xpetra::EpetraCrsGraphT\n"
466 "Use node_t that is supported by "
467 "Xpetra with Epetra classes\n")
468 + e.what());
469 }
470 return xa;
471 }
472
473 static RCP<Epetra_CrsGraph> doMigration(const Epetra_CrsGraph &from,
474 size_t numLocalRows, const gno_t *myNewRows)
475 {
476 // source map
477 const Epetra_BlockMap &smap = from.RowMap();
478 gno_t numGlobalRows = smap.NumGlobalElements();
479 lno_t oldNumElts = smap.NumMyElements();
480 int base = smap.MinAllGID();
481
482 // target map
483 const Epetra_Comm &comm = from.Comm();
484 Epetra_BlockMap tmap(numGlobalRows, numLocalRows, myNewRows, 1, base, comm);
485 lno_t newNumElts = tmap.NumMyElements();
486
487 // importer
488 Epetra_Import importer(tmap, smap);
489
490 // number of non zeros in my new rows
491 Epetra_Vector numOld(smap);
492 Epetra_Vector numNew(tmap);
493
494 for (int lid=0; lid < oldNumElts; lid++){
495 numOld[lid] = from.NumMyIndices(lid);
496 }
497 numNew.Import(numOld, importer, Insert);
498
499 Array<int> nnz(newNumElts);
500 for (int lid=0; lid < newNumElts; lid++){
501 nnz[lid] = static_cast<int>(numNew[lid]);
502 }
503
504 // target graph
505 RCP<Epetra_CrsGraph> G = rcp(new Epetra_CrsGraph(::Copy, tmap, nnz.getRawPtr(), true));
506 G->Import(from, importer, Insert);
507 G->FillComplete();
508
509 return G;
510 }
511
512};
513#endif
514
516// Xpetra::CrsGraph
517// KDDKDD Do we need specializations for Xpetra::TpetraCrsGraph and
518// KDDKDD Xpetra::EpetraCrsGraph?
519template <typename lno_t,
520 typename gno_t,
521 typename node_t>
522struct XpetraTraits<Xpetra::CrsGraph<lno_t, gno_t, node_t> >
523{
524 typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
525 typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
526 typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
527
528 static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
529 {
530 return a;
531 }
532
533 static RCP<x_graph_t> doMigration(const x_graph_t &from,
534 size_t numLocalRows, const gno_t *myNewRows)
535 {
536 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
537
538 if (lib == Xpetra::UseEpetra){
539 throw std::logic_error("compiler should have used specialization");
540 } else{
541 // Do the import with the Tpetra::CrsGraph traits object
542 const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
543 RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
544
545 RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
546 *tg, numLocalRows, myNewRows);
547
548 RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
549 return xgnew;
550 }
551 }
552};
553
554
556// Xpetra::CrsGraph specialization
557template < typename node_t>
558struct XpetraTraits<Xpetra::CrsGraph<int, int, node_t> >
559{
560 typedef int lno_t;
561 typedef int gno_t;
562 typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
563 typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
564 typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
565
566 static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
567 {
568 return a;
569 }
570
571 static RCP<x_graph_t> doMigration(const x_graph_t &from,
572 size_t numLocalRows, const gno_t *myNewRows)
573 {
574 Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
575
576 if (lib == Xpetra::UseEpetra){
577#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
578 typedef Xpetra::EpetraCrsGraphT<gno_t,node_t> xe_graph_t;
579 typedef Epetra_CrsGraph e_graph_t;
580 // Do the import with the Epetra_CrsGraph traits object
581 const xe_graph_t *xeg = dynamic_cast<const xe_graph_t *>(&from);
582 RCP<const e_graph_t> eg = xeg->getEpetra_CrsGraph();
583
584 RCP<e_graph_t> egnew = XpetraTraits<e_graph_t>::doMigration(
585 *eg, numLocalRows, myNewRows);
586
587 RCP<x_graph_t> xgnew = XpetraTraits<e_graph_t>::convertToXpetra(egnew);
588
589 return xgnew;
590#else
591 throw std::runtime_error("Xpetra with Epetra requested, but "
592 "Trilinos is not built with Epetra");
593#endif
594 } else{
595 // Do the import with the Tpetra::CrsGraph traits object
596 const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
597 RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
598
599 RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
600 *tg, numLocalRows, myNewRows);
601
602 RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
603
604 return xgnew;
605 }
606 }
607};
608
610// Tpetra::Vector
611template <typename scalar_t,
612 typename lno_t,
613 typename gno_t,
614 typename node_t>
615struct XpetraTraits<Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
616{
617 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
618 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
619 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
620
621 static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
622 {
623 return rcp(new xt_vector_t(a));
624 }
625
626 static RCP<t_vector_t> doMigration(const t_vector_t &from,
627 size_t numLocalElts, const gno_t *myNewElts)
628 {
629 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
630
631 // source map
632 const RCP<const map_t> &smap = from.getMap();
633 gno_t numGlobalElts = smap->getGlobalNumElements();
634 gno_t base = smap->getMinAllGlobalIndex();
635
636 // target map
637 ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
638 const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
639 RCP<const map_t> tmap = rcp(new map_t(numGlobalElts, eltList, base, comm));
640
641 // importer
642 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
643
644 // target vector
645 RCP<t_vector_t> V =
646 Tpetra::createVector<scalar_t,lno_t,gno_t,node_t>(tmap);
647 V->doImport(from, importer, Tpetra::INSERT);
648
649 return V;
650 }
651};
652
654#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
655// Epetra_Vector
656template < >
657struct XpetraTraits<Epetra_Vector>
658{
663
664 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
665
666 static inline RCP<x_vector_t> convertToXpetra(const RCP<Epetra_Vector> &a)
667 {
668 RCP<Xpetra::EpetraVectorT<gno_t, node_t> > xev;
669 try {
670 xev = rcp(new Xpetra::EpetraVectorT<gno_t,node_t>(a));
671 }
672 catch (std::exception &e) {
673 if (std::is_same<node_t, Xpetra::EpetraNode>::value)
674 throw std::runtime_error(std::string("Cannot convert from "
675 "Epetra_Vector to "
676 "Xpetra::EpetraVectorT\n")
677 + e.what());
678 else
679 throw std::runtime_error(std::string("Cannot convert from "
680 "Epetra_Vector to "
681 "Xpetra::EpetraVectorT\n"
682 "Use node_t that is supported by "
683 "Xpetra with Epetra classes\n")
684 + e.what());
685 }
686 return rcp_implicit_cast<x_vector_t>(xev);
687 }
688
689 static RCP<Epetra_Vector> doMigration(const Epetra_Vector &from,
690 size_t numLocalElts, const gno_t *myNewElts)
691 {
692 // source map
693 const Epetra_BlockMap &smap = from.Map();
694 gno_t numGlobalElts = smap.NumGlobalElements();
695 int base = smap.MinAllGID();
696
697 // target map
698 const Epetra_Comm &comm = from.Comm();
699 const Epetra_BlockMap tmap(numGlobalElts, numLocalElts, myNewElts,
700 1, base, comm);
701
702 // importer
703 Epetra_Import importer(tmap, smap);
704
705 // target vector
706 RCP<Epetra_Vector> V = rcp(new Epetra_Vector(tmap, true));
707 Epetra_CombineMode c = Insert;
708 V->Import(from, importer, c);
709
710 return V;
711 }
712};
713#endif
714
716// Xpetra::Vector
717template <typename scalar_t,
718 typename lno_t,
719 typename gno_t,
720 typename node_t>
721struct XpetraTraits<Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
722{
723 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
724 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
725 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
726
727 static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
728 {
729 return a;
730 }
731
732 static RCP<x_vector_t> doMigration(const x_vector_t &from,
733 size_t numLocalRows, const gno_t *myNewRows)
734 {
735 Xpetra::UnderlyingLib lib = from.getMap()->lib();
736
737 if (lib == Xpetra::UseEpetra){
738 throw std::logic_error("compiler should have used specialization");
739 } else{
740 // Do the import with the Tpetra::Vector traits object
741 const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
742 RCP<const t_vector_t> tv = xtv->getTpetra_Vector();
743
744 RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
745 *tv, numLocalRows, myNewRows);
746
747 RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
748
749 return xvnew;
750 }
751 }
752};
753
755// Xpetra::Vector specialization
756template <typename node_t>
757struct XpetraTraits<Xpetra::Vector<double, int, int, node_t> >
758{
759 typedef double scalar_t;
760 typedef int lno_t;
761 typedef int gno_t;
762 typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
763 typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
764 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
765
766 static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
767 {
768 return a;
769 }
770
771 static RCP<x_vector_t> doMigration(const x_vector_t &from,
772 size_t numLocalRows, const gno_t *myNewRows)
773 {
774 Xpetra::UnderlyingLib lib = from.getMap()->lib();
775
776 if (lib == Xpetra::UseEpetra){
777#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
778 typedef Epetra_Vector e_vector_t;
779 typedef Xpetra::EpetraVectorT<gno_t,node_t> xe_vector_t;
780 // Do the import with the Epetra_Vector traits object
781 const xe_vector_t *xev = dynamic_cast<const xe_vector_t *>(&from);
782 RCP<const e_vector_t> ev = rcp(xev->getEpetra_Vector());
783
784 RCP<e_vector_t> evnew = XpetraTraits<e_vector_t>::doMigration(
785 *ev, numLocalRows, myNewRows);
786
787 RCP<x_vector_t> xvnew = XpetraTraits<e_vector_t>::convertToXpetra(evnew);
788
789 return xvnew;
790#else
791 throw std::runtime_error("Xpetra with Epetra requested, but "
792 "Trilinos is not built with Epetra");
793#endif
794 } else{
795 // Do the import with the Tpetra::Vector traits object
796 const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
797 RCP<t_vector_t> tv = xtv->getTpetra_Vector();
798
799 RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
800 *tv, numLocalRows, myNewRows);
801
802 RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
803
804 return xvnew;
805 }
806 }
807};
808
810// Tpetra::MultiVector
811template <typename scalar_t,
812 typename lno_t,
813 typename gno_t,
814 typename node_t>
815struct XpetraTraits<Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
816{
817 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
818 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
819 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
820
821 static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
822 {
823 return rcp(new xt_vector_t(a));
824 }
825
826 static RCP<t_vector_t> doMigration(const t_vector_t &from,
827 size_t numLocalElts, const gno_t *myNewElts)
828 {
829 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
830
831 // source map
832 const RCP<const map_t> &smap = from.getMap();
833 gno_t numGlobalElts = smap->getGlobalNumElements();
834 gno_t base = smap->getMinAllGlobalIndex();
835
836 // target map
837 ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
838 const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
839 RCP<const map_t> tmap = rcp(new map_t(numGlobalElts, eltList, base, comm));
840
841 // importer
842 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
843
844 // target vector
845 RCP<t_vector_t> MV = rcp(
846 new t_vector_t(tmap, from.getNumVectors(), true));
847 MV->doImport(from, importer, Tpetra::INSERT);
848
849 return MV;
850 }
851};
852
854#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
855// Epetra_MultiVector
856template < >
857struct XpetraTraits<Epetra_MultiVector>
858{
863 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
864
865 static inline RCP<x_mvector_t> convertToXpetra(
866 const RCP<Epetra_MultiVector> &a)
867 {
868 RCP<Xpetra::EpetraMultiVectorT<gno_t, node_t> > xemv;
869 try {
870 xemv = rcp(new Xpetra::EpetraMultiVectorT<gno_t,node_t>(a));
871 }
872 catch (std::exception &e) {
873 if (std::is_same<node_t, Xpetra::EpetraNode>::value)
874 throw std::runtime_error(std::string("Cannot convert from "
875 "Epetra_MultiVector to "
876 "Xpetra::EpetraMultiVectorT\n")
877 + e.what());
878 else
879 throw std::runtime_error(std::string("Cannot convert from "
880 "Epetra_MultiVector to "
881 "Xpetra::EpetraMultiVectorT\n"
882 "Use node_t that is supported by "
883 "Xpetra with Epetra classes\n")
884 + e.what());
885 }
886 return rcp_implicit_cast<x_mvector_t>(xemv);
887 }
888
889 static RCP<Epetra_MultiVector> doMigration(const Epetra_MultiVector &from,
890 size_t numLocalElts, const gno_t *myNewElts)
891 {
892 // source map
893 const Epetra_BlockMap &smap = from.Map();
894 gno_t numGlobalElts = smap.NumGlobalElements();
895 int base = smap.MinAllGID();
896
897 // target map
898 const Epetra_Comm &comm = from.Comm();
899 const Epetra_BlockMap tmap(numGlobalElts, numLocalElts, myNewElts,
900 1, base, comm);
901
902 // importer
903 Epetra_Import importer(tmap, smap);
904
905 // target vector
906 RCP<Epetra_MultiVector> MV = rcp(
907 new Epetra_MultiVector(tmap, from.NumVectors(), true));
908 Epetra_CombineMode c = Insert;
909 MV->Import(from, importer, c);
910
911 return MV;
912 }
913};
914#endif
915
917// Xpetra::MultiVector
918template <typename scalar_t,
919 typename lno_t,
920 typename gno_t,
921 typename node_t>
922struct XpetraTraits<Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
923{
924 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
925 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
926 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
927
928 static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
929 {
930 return a;
931 }
932
933 static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
934 size_t numLocalRows, const gno_t *myNewRows)
935 {
936 Xpetra::UnderlyingLib lib = from.getMap()->lib();
937
938 if (lib == Xpetra::UseEpetra){
939 throw std::logic_error("compiler should have used specialization");
940 } else{
941 // Do the import with the Tpetra::MultiVector traits object
942 const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
943 RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
944
945 RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
946 *tv, numLocalRows, myNewRows);
947
948 RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
949
950 return xvnew;
951 }
952 }
953};
954
956// Xpetra::MultiVector specialization
957template <typename node_t>
958struct XpetraTraits<Xpetra::MultiVector<double, int, int, node_t> >
959{
960 typedef double scalar_t;
961 typedef int lno_t;
962 typedef int gno_t;
963 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
964 typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
965 typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
966
967 static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
968 {
969 return a;
970 }
971
972 static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
973 size_t numLocalRows, const gno_t *myNewRows)
974 {
975 Xpetra::UnderlyingLib lib = from.getMap()->lib();
976
977 if (lib == Xpetra::UseEpetra){
978#if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
979 typedef Epetra_MultiVector e_mvector_t;
980 typedef Xpetra::EpetraMultiVectorT<gno_t,node_t> xe_mvector_t;
981 // Do the import with the Epetra_MultiVector traits object
982 const xe_mvector_t *xev = dynamic_cast<const xe_mvector_t *>(&from);
983 RCP<e_mvector_t> ev = xev->getEpetra_MultiVector();
984
985 RCP<e_mvector_t> evnew = XpetraTraits<e_mvector_t>::doMigration(
986 *ev, numLocalRows, myNewRows);
987
988 RCP<x_mvector_t> xvnew = XpetraTraits<e_mvector_t>::convertToXpetra(evnew);
989
990 return xvnew;
991#else
992 throw std::runtime_error("Xpetra with Epetra requested, but "
993 "Trilinos is not built with Epetra");
994#endif
995 } else{
996 // Do the import with the Tpetra::MultiVector traits object
997 const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
998 RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
999
1000 RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
1001 *tv, numLocalRows, myNewRows);
1002
1003 RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
1004
1005 return xvnew;
1006 }
1007 }
1008};
1009
1010#endif // DOXYGEN_SHOULD_SKIP_THIS
1011
1012} //namespace Zoltan2
1013
1014#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
default_gno_t gno_t
The ordinal type (e.g., int, long, int64_t) that can represent global counts and identifiers.
default_node_t node_t
The Kokkos node type. This is only meaningful for users of Tpetra objects.
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices.
default_scalar_t scalar_t
The data type for weights and coordinates.
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.