Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_XpetraTraits.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
50 #ifndef _ZOLTAN2_XPETRATRAITS_HPP_
51 #define _ZOLTAN2_XPETRATRAITS_HPP_
52 
53 #include <Zoltan2_InputTraits.hpp>
54 #include <Zoltan2_Standards.hpp>
55 
56 #include <Xpetra_TpetraCrsMatrix.hpp>
57 #include <Xpetra_TpetraVector.hpp>
58 #include <Tpetra_Vector.hpp>
59 
60 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
61 #include <Xpetra_EpetraCrsMatrix.hpp>
62 #include <Xpetra_EpetraVector.hpp>
63 #include <Xpetra_EpetraUtils.hpp>
64 #endif
65 
66 namespace Zoltan2 {
67 
69 // Extra traits needed only for Xpetra matrices and graphs
70 
93 template <typename User>
95 {
98  static inline RCP<User> convertToXpetra(const RCP<User> &a)
99  {
100  return a;
101  }
102 
106 
110 
116  static RCP<User> doMigration(const User &from,
117  size_t numLocalRows, const gno_t *myNewRows)
118  {
119  return Teuchos::null;
120  }
121 };
122 
123 #ifndef DOXYGEN_SHOULD_SKIP_THIS
124 
126 // Tpetra::CrsMatrix
127 template <typename scalar_t,
128  typename lno_t,
129  typename gno_t,
130  typename node_t>
131 struct XpetraTraits<Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
132 {
133  typedef typename Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> xmatrix_t;
134  typedef typename Xpetra::TpetraCrsMatrix<scalar_t,lno_t,gno_t,node_t> xtmatrix_t;
135  typedef typename Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> tmatrix_t;
136 
137  static inline RCP<xmatrix_t> convertToXpetra(const RCP<tmatrix_t> &a)
138  {
139  return rcp(new xtmatrix_t(a));
140  }
141 
142  static RCP<tmatrix_t> doMigration(const tmatrix_t &from,
143  size_t numLocalRows, const gno_t *myNewRows)
144  {
145  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
146 
147  // source map
148  const RCP<const map_t> &smap = from.getRowMap();
149  gno_t numGlobalRows = smap->getGlobalNumElements();
150  gno_t base = smap->getMinAllGlobalIndex();
151 
152  // target map
153  ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
154  const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
155  RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
156 
157  // importer
158  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
159 
160  // target matrix
161  // Chris Siefert proposed using the following to make migration
162  // more efficient.
163  // By default, the Domain and Range maps are the same as in "from".
164  // As in the original code, we instead set them both to tmap.
165  // The assumption is a square matrix.
166  // TODO: what about rectangular matrices?
167  // TODO: Should choice of domain/range maps be an option to this function?
168 
169  // KDD 3/7/16: disabling Chris' new code to avoid dashboard failures;
170  // KDD 3/7/16: can re-enable when issue #114 is fixed.
171  // KDD 3/7/16: when re-enable CSIEFERT code, can comment out
172  // KDD 3/7/16: "Original way" code.
173  // KDD 1/27/17: Re-enabling Chris' code, as this issue is resolved.
174  RCP<tmatrix_t> M;
175  from.importAndFillComplete(M, importer, tmap, tmap);
176 
179  //int oldNumElts = smap->getNodeNumElements();
180  //int newNumElts = numLocalRows;
181 
183  //typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
184  //vector_t numOld(smap); // TODO These vectors should have scalar=size_t,
185  //vector_t numNew(tmap); // but ETI does not yet support that.
186  //for (int lid=0; lid < oldNumElts; lid++){
187  //numOld.replaceGlobalValue(smap->getGlobalElement(lid),
188  //scalar_t(from.getNumEntriesInLocalRow(lid)));
189  //}
190  //numNew.doImport(numOld, importer, Tpetra::INSERT);
191 
192  // TODO Could skip this copy if could declare vector with scalar=size_t.
193  //ArrayRCP<size_t> nnz(newNumElts);
194  //if (newNumElts > 0){
195  //ArrayRCP<scalar_t> ptr = numNew.getDataNonConst(0);
196  //for (int lid=0; lid < newNumElts; lid++){
197  //nnz[lid] = static_cast<size_t>(ptr[lid]);
198  //}
199  //}
200 
201  //RCP<tmatrix_t> M = rcp(new tmatrix_t(tmap, nnz, Tpetra::StaticProfile));
202  //M->doImport(from, importer, Tpetra::INSERT);
203  //M->fillComplete();
204 
205  // End of original way we did it.
206  return M;
207  }
208 };
209 
211 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
212 // Epetra_CrsMatrix
213 template <>
214 struct XpetraTraits<Epetra_CrsMatrix>
215 {
216  typedef InputTraits<Epetra_CrsMatrix>::scalar_t scalar_t;
217  typedef InputTraits<Epetra_CrsMatrix>::lno_t lno_t;
218  typedef InputTraits<Epetra_CrsMatrix>::gno_t gno_t;
219  typedef InputTraits<Epetra_CrsMatrix>::node_t node_t;
220 
221  static inline RCP<Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> >
222  convertToXpetra(const RCP<Epetra_CrsMatrix> &a)
223  {
224  RCP<Xpetra::EpetraCrsMatrixT<gno_t, node_t> > xa;
225  try {
226  xa = rcp(new Xpetra::EpetraCrsMatrixT<gno_t, node_t>(a));
227  }
228  catch (std::exception &e) {
229  if (std::is_same<node_t, Xpetra::EpetraNode>::value)
230  throw std::runtime_error(std::string("Cannot convert from "
231  "Epetra_CrsMatrix to "
232  "Xpetra::EpetraCrsMatrixT\n")
233  + e.what());
234  else
235  throw std::runtime_error(std::string("Cannot convert from "
236  "Epetra_CrsMatrix to "
237  "Xpetra::EpetraCrsMatrixT\n"
238  "Use node_t that is supported by "
239  "Xpetra with Epetra classes\n")
240  + e.what());
241  }
242  return xa;
243  }
244 
245 
246  static RCP<Epetra_CrsMatrix> doMigration(const Epetra_CrsMatrix &from,
247  size_t numLocalRows, const gno_t *myNewRows)
248  {
249  // source map
250  const Epetra_Map &smap = from.RowMap();
251  gno_t numGlobalRows = smap.NumGlobalElements();
252  int base = smap.MinAllGID();
253 
254  // target map
255  const Epetra_Comm &comm = from.Comm();
256  Epetra_Map tmap(numGlobalRows, numLocalRows, myNewRows, base, comm);
257 
258  // importer
259  Epetra_Import importer(tmap, smap);
260 
261 
262  // target matrix
263  // Chris Siefert proposed using the following to make migration
264  // more efficient.
265  // By default, the Domain and Range maps are the same as in "from".
266  // As in the original code, we instead set them both to tmap.
267  // The assumption is a square matrix.
268  // TODO: what about rectangular matrices?
269  // TODO: Should choice of domain/range maps be an option to this function?
270 
271  RCP<Epetra_CrsMatrix> M = rcp(new Epetra_CrsMatrix(from, importer,
272  &tmap, &tmap));
273 
274  // Original way we did it:
275  //
276  // int oldNumElts = smap.NumMyElements();
277  //
278  // // number of non zeros in my new rows
279  // Epetra_Vector numOld(smap);
280  // Epetra_Vector numNew(tmap);
281  //
282  // for (int lid=0; lid < oldNumElts; lid++){
283  // numOld[lid] = from.NumMyEntries(lid);
284  // }
285  // numNew.Import(numOld, importer, Insert);
286  //
287  // int newNumElts = numLocalRows;
288  // Array<int> nnz(newNumElts);
289  // for (int lid=0; lid < newNumElts; lid++){
290  // nnz[lid] = static_cast<int>(numNew[lid]);
291  // }
292  //
293  // RCP<Epetra_CrsMatrix> M =
294  // rcp(new Epetra_CrsMatrix(::Copy, tmap, nnz.getRawPtr(), true));
295  // M->Import(from, importer, Insert);
296  // M->FillComplete();
297 
298  return M;
299  }
300 };
301 #endif
302 
304 // Xpetra::CrsMatrix
305 // KDDKDD: Do we need specializations for Xpetra::EpetraCrsMatrix and
306 // KDDKDD: Xpetra::TpetraCrsMatrix
307 template <typename scalar_t,
308  typename lno_t,
309  typename gno_t,
310  typename node_t>
311 struct XpetraTraits<Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
312 {
313  typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
314  typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
315  typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
316 
317  static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
318  {
319  return a;
320  }
321 
322  static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
323  size_t numLocalRows, const gno_t *myNewRows)
324  {
325  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
326 
327  if (lib == Xpetra::UseEpetra){
328  throw std::logic_error("compiler should have used specialization");
329  } else{
330  // Do the import with the Tpetra::CrsMatrix traits object
331  const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
332  RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
333 
334  RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
335  *tm, numLocalRows, myNewRows);
336 
337  RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
338 
339  return xmnew;
340  }
341  }
342 };
343 
345 // Xpetra::CrsMatrix specialization
346 
347 template <typename node_t>
348 struct XpetraTraits<Xpetra::CrsMatrix<double, int, int, node_t> >
349 {
350  typedef double scalar_t;
351  typedef int lno_t;
352  typedef int gno_t;
353  typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
354  typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
355  typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
356 
357  static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
358  {
359  return a;
360  }
361 
362  static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
363  size_t numLocalRows, const gno_t *myNewRows)
364  {
365  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
366 
367  if (lib == Xpetra::UseEpetra){
368 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
369  typedef Epetra_CrsMatrix e_matrix_t;
370  typedef Xpetra::EpetraCrsMatrixT<gno_t,node_t> xe_matrix_t;
371  // Do the import with the Epetra_CrsMatrix traits object
372  const xe_matrix_t *xem = dynamic_cast<const xe_matrix_t *>(&from);
373  RCP<const e_matrix_t> em = xem->getEpetra_CrsMatrix();
374 
375  RCP<e_matrix_t> emnew = XpetraTraits<e_matrix_t>::doMigration(
376  *em, numLocalRows, myNewRows);
377 
378  RCP<x_matrix_t> xmnew = XpetraTraits<e_matrix_t>::convertToXpetra(emnew);
379 
380  return xmnew;
381 #else
382  throw std::runtime_error("Xpetra with Epetra requested, but "
383  "Trilinos is not built with Epetra");
384 #endif
385  } else{
386  // Do the import with the Tpetra::CrsMatrix traits object
387  const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
388  RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
389 
390  RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
391  *tm, numLocalRows, myNewRows);
392 
393  RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
394 
395  return xmnew;
396  }
397  }
398 };
399 
400 
402 // Tpetra::CrsGraph
403 template <typename lno_t,
404  typename gno_t,
405  typename node_t>
406 struct XpetraTraits<Tpetra::CrsGraph<lno_t, gno_t, node_t> >
407 {
408  typedef typename Xpetra::CrsGraph<lno_t, gno_t, node_t> xgraph_t;
409  typedef typename Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xtgraph_t;
410  typedef typename Tpetra::CrsGraph<lno_t, gno_t, node_t> tgraph_t;
411 
412  static inline RCP<xgraph_t> convertToXpetra(const RCP<tgraph_t> &a)
413  {
414  return rcp(new xtgraph_t(a));
415  }
416 
417  static RCP<tgraph_t> doMigration(const tgraph_t &from,
418  size_t numLocalRows, const gno_t *myNewRows)
419  {
420  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
421 
422  // source map
423  const RCP<const map_t> &smap = from.getRowMap();
424  int oldNumElts = smap->getNodeNumElements();
425  gno_t numGlobalRows = smap->getGlobalNumElements();
426  gno_t base = smap->getMinAllGlobalIndex();
427 
428  // target map
429  ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
430  const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
431  RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
432 
433  // importer
434  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
435 
436  // number of entries in my new rows
437  typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
438  vector_t numOld(smap);
439  vector_t numNew(tmap);
440  for (int lid=0; lid < oldNumElts; lid++){
441  numOld.replaceGlobalValue(smap->getGlobalElement(lid),
442  from.getNumEntriesInLocalRow(lid));
443  }
444  numNew.doImport(numOld, importer, Tpetra::INSERT);
445 
446  size_t numElts = tmap->getNodeNumElements();
447  ArrayRCP<const gno_t> nnz;
448  if (numElts > 0)
449  nnz = numNew.getData(0); // hangs if vector len == 0
450 
451  ArrayRCP<const size_t> nnz_size_t;
452 
453  if (numElts && sizeof(gno_t) != sizeof(size_t)){
454  size_t *vals = new size_t [numElts];
455  nnz_size_t = arcp(vals, 0, numElts, true);
456  for (size_t i=0; i < numElts; i++){
457  vals[i] = static_cast<size_t>(nnz[i]);
458  }
459  }
460  else{
461  nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
462  }
463 
464  // target graph
465  RCP<tgraph_t> G = rcp(new tgraph_t(tmap, nnz_size_t(), Tpetra::StaticProfile));
466 
467  G->doImport(from, importer, Tpetra::INSERT);
468  G->fillComplete();
469 
470  return G;
471  }
472 
473 };
474 
476 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
477 // Epetra_CrsGraph
478 template < >
479 struct XpetraTraits<Epetra_CrsGraph>
480 {
481  typedef InputTraits<Epetra_CrsGraph>::lno_t lno_t;
482  typedef InputTraits<Epetra_CrsGraph>::gno_t gno_t;
483  typedef InputTraits<Epetra_CrsGraph>::node_t node_t;
484  static inline RCP<Xpetra::CrsGraph<lno_t,gno_t,node_t> >
485  convertToXpetra(const RCP<Epetra_CrsGraph> &a)
486  {
487  RCP<Xpetra::EpetraCrsGraphT<gno_t, node_t> > xa;
488  try {
489  xa = rcp(new Xpetra::EpetraCrsGraphT<gno_t, node_t>(a));
490  }
491  catch (std::exception &e) {
492  if (std::is_same<node_t, Xpetra::EpetraNode>::value)
493  throw std::runtime_error(std::string("Cannot convert from "
494  "Epetra_CrsGraph to "
495  "Xpetra::EpetraCrsGraphT\n")
496  + e.what());
497  else
498  throw std::runtime_error(std::string("Cannot convert from "
499  "Epetra_CrsGraph to "
500  "Xpetra::EpetraCrsGraphT\n"
501  "Use node_t that is supported by "
502  "Xpetra with Epetra classes\n")
503  + e.what());
504  }
505  return xa;
506  }
507 
508  static RCP<Epetra_CrsGraph> doMigration(const Epetra_CrsGraph &from,
509  size_t numLocalRows, const gno_t *myNewRows)
510  {
511  // source map
512  const Epetra_BlockMap &smap = from.RowMap();
513  gno_t numGlobalRows = smap.NumGlobalElements();
514  lno_t oldNumElts = smap.NumMyElements();
515  int base = smap.MinAllGID();
516 
517  // target map
518  const Epetra_Comm &comm = from.Comm();
519  Epetra_BlockMap tmap(numGlobalRows, numLocalRows, myNewRows, 1, base, comm);
520  lno_t newNumElts = tmap.NumMyElements();
521 
522  // importer
523  Epetra_Import importer(tmap, smap);
524 
525  // number of non zeros in my new rows
526  Epetra_Vector numOld(smap);
527  Epetra_Vector numNew(tmap);
528 
529  for (int lid=0; lid < oldNumElts; lid++){
530  numOld[lid] = from.NumMyIndices(lid);
531  }
532  numNew.Import(numOld, importer, Insert);
533 
534  Array<int> nnz(newNumElts);
535  for (int lid=0; lid < newNumElts; lid++){
536  nnz[lid] = static_cast<int>(numNew[lid]);
537  }
538 
539  // target graph
540  RCP<Epetra_CrsGraph> G = rcp(new Epetra_CrsGraph(::Copy, tmap, nnz.getRawPtr(), true));
541  G->Import(from, importer, Insert);
542  G->FillComplete();
543 
544  return G;
545  }
546 
547 };
548 #endif
549 
551 // Xpetra::CrsGraph
552 // KDDKDD Do we need specializations for Xpetra::TpetraCrsGraph and
553 // KDDKDD Xpetra::EpetraCrsGraph?
554 template <typename lno_t,
555  typename gno_t,
556  typename node_t>
557 struct XpetraTraits<Xpetra::CrsGraph<lno_t, gno_t, node_t> >
558 {
559  typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
560  typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
561  typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
562 
563  static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
564  {
565  return a;
566  }
567 
568  static RCP<x_graph_t> doMigration(const x_graph_t &from,
569  size_t numLocalRows, const gno_t *myNewRows)
570  {
571  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
572 
573  if (lib == Xpetra::UseEpetra){
574  throw std::logic_error("compiler should have used specialization");
575  } else{
576  // Do the import with the Tpetra::CrsGraph traits object
577  const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
578  RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
579 
580  RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
581  *tg, numLocalRows, myNewRows);
582 
583  RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
584  return xgnew;
585  }
586  }
587 };
588 
589 
591 // Xpetra::CrsGraph specialization
592 template < typename node_t>
593 struct XpetraTraits<Xpetra::CrsGraph<int, int, node_t> >
594 {
595  typedef int lno_t;
596  typedef int gno_t;
597  typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
598  typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
599  typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
600 
601  static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
602  {
603  return a;
604  }
605 
606  static RCP<x_graph_t> doMigration(const x_graph_t &from,
607  size_t numLocalRows, const gno_t *myNewRows)
608  {
609  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
610 
611  if (lib == Xpetra::UseEpetra){
612 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
613  typedef Xpetra::EpetraCrsGraphT<gno_t,node_t> xe_graph_t;
614  typedef Epetra_CrsGraph e_graph_t;
615  // Do the import with the Epetra_CrsGraph traits object
616  const xe_graph_t *xeg = dynamic_cast<const xe_graph_t *>(&from);
617  RCP<const e_graph_t> eg = xeg->getEpetra_CrsGraph();
618 
619  RCP<e_graph_t> egnew = XpetraTraits<e_graph_t>::doMigration(
620  *eg, numLocalRows, myNewRows);
621 
622  RCP<x_graph_t> xgnew = XpetraTraits<e_graph_t>::convertToXpetra(egnew);
623 
624  return xgnew;
625 #else
626  throw std::runtime_error("Xpetra with Epetra requested, but "
627  "Trilinos is not built with Epetra");
628 #endif
629  } else{
630  // Do the import with the Tpetra::CrsGraph traits object
631  const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
632  RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
633 
634  RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
635  *tg, numLocalRows, myNewRows);
636 
637  RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
638 
639  return xgnew;
640  }
641  }
642 };
643 
645 // Tpetra::Vector
646 template <typename scalar_t,
647  typename lno_t,
648  typename gno_t,
649  typename node_t>
650 struct XpetraTraits<Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
651 {
652  typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
653  typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
654  typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
655 
656  static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
657  {
658  return rcp(new xt_vector_t(a));
659  }
660 
661  static RCP<t_vector_t> doMigration(const t_vector_t &from,
662  size_t numLocalElts, const gno_t *myNewElts)
663  {
664  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
665 
666  // source map
667  const RCP<const map_t> &smap = from.getMap();
668  gno_t numGlobalElts = smap->getGlobalNumElements();
669  gno_t base = smap->getMinAllGlobalIndex();
670 
671  // target map
672  ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
673  const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
674  RCP<const map_t> tmap = rcp(new map_t(numGlobalElts, eltList, base, comm));
675 
676  // importer
677  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
678 
679  // target vector
680  RCP<t_vector_t> V =
681  Tpetra::createVector<scalar_t,lno_t,gno_t,node_t>(tmap);
682  V->doImport(from, importer, Tpetra::INSERT);
683 
684  return V;
685  }
686 };
687 
689 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
690 // Epetra_Vector
691 template < >
692 struct XpetraTraits<Epetra_Vector>
693 {
694  typedef InputTraits<Epetra_Vector>::lno_t lno_t;
695  typedef InputTraits<Epetra_Vector>::gno_t gno_t;
696  typedef InputTraits<Epetra_Vector>::node_t node_t;
697  typedef InputTraits<Epetra_Vector>::scalar_t scalar_t;
698 
699  typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
700 
701  static inline RCP<x_vector_t> convertToXpetra(const RCP<Epetra_Vector> &a)
702  {
703  RCP<Xpetra::EpetraVectorT<gno_t, node_t> > xev;
704  try {
705  xev = rcp(new Xpetra::EpetraVectorT<gno_t,node_t>(a));
706  }
707  catch (std::exception &e) {
708  if (std::is_same<node_t, Xpetra::EpetraNode>::value)
709  throw std::runtime_error(std::string("Cannot convert from "
710  "Epetra_Vector to "
711  "Xpetra::EpetraVectorT\n")
712  + e.what());
713  else
714  throw std::runtime_error(std::string("Cannot convert from "
715  "Epetra_Vector to "
716  "Xpetra::EpetraVectorT\n"
717  "Use node_t that is supported by "
718  "Xpetra with Epetra classes\n")
719  + e.what());
720  }
721  return rcp_implicit_cast<x_vector_t>(xev);
722  }
723 
724  static RCP<Epetra_Vector> doMigration(const Epetra_Vector &from,
725  size_t numLocalElts, const gno_t *myNewElts)
726  {
727  // source map
728  const Epetra_BlockMap &smap = from.Map();
729  gno_t numGlobalElts = smap.NumGlobalElements();
730  int base = smap.MinAllGID();
731 
732  // target map
733  const Epetra_Comm &comm = from.Comm();
734  const Epetra_BlockMap tmap(numGlobalElts, numLocalElts, myNewElts,
735  1, base, comm);
736 
737  // importer
738  Epetra_Import importer(tmap, smap);
739 
740  // target vector
741  RCP<Epetra_Vector> V = rcp(new Epetra_Vector(tmap, true));
742  Epetra_CombineMode c = Insert;
743  V->Import(from, importer, c);
744 
745  return V;
746  }
747 };
748 #endif
749 
751 // Xpetra::Vector
752 template <typename scalar_t,
753  typename lno_t,
754  typename gno_t,
755  typename node_t>
756 struct XpetraTraits<Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
757 {
758  typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
759  typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
760  typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
761 
762  static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
763  {
764  return a;
765  }
766 
767  static RCP<x_vector_t> doMigration(const x_vector_t &from,
768  size_t numLocalRows, const gno_t *myNewRows)
769  {
770  Xpetra::UnderlyingLib lib = from.getMap()->lib();
771 
772  if (lib == Xpetra::UseEpetra){
773  throw std::logic_error("compiler should have used specialization");
774  } else{
775  // Do the import with the Tpetra::Vector traits object
776  const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
777  RCP<const t_vector_t> tv = xtv->getTpetra_Vector();
778 
779  RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
780  *tv, numLocalRows, myNewRows);
781 
782  RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
783 
784  return xvnew;
785  }
786  }
787 };
788 
790 // Xpetra::Vector specialization
791 template <typename node_t>
792 struct XpetraTraits<Xpetra::Vector<double, int, int, node_t> >
793 {
794  typedef double scalar_t;
795  typedef int lno_t;
796  typedef int gno_t;
797  typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
798  typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
799  typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
800 
801  static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
802  {
803  return a;
804  }
805 
806  static RCP<x_vector_t> doMigration(const x_vector_t &from,
807  size_t numLocalRows, const gno_t *myNewRows)
808  {
809  Xpetra::UnderlyingLib lib = from.getMap()->lib();
810 
811  if (lib == Xpetra::UseEpetra){
812 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
813  typedef Epetra_Vector e_vector_t;
814  typedef Xpetra::EpetraVectorT<gno_t,node_t> xe_vector_t;
815  // Do the import with the Epetra_Vector traits object
816  const xe_vector_t *xev = dynamic_cast<const xe_vector_t *>(&from);
817  RCP<const e_vector_t> ev = rcp(xev->getEpetra_Vector());
818 
819  RCP<e_vector_t> evnew = XpetraTraits<e_vector_t>::doMigration(
820  *ev, numLocalRows, myNewRows);
821 
822  RCP<x_vector_t> xvnew = XpetraTraits<e_vector_t>::convertToXpetra(evnew);
823 
824  return xvnew;
825 #else
826  throw std::runtime_error("Xpetra with Epetra requested, but "
827  "Trilinos is not built with Epetra");
828 #endif
829  } else{
830  // Do the import with the Tpetra::Vector traits object
831  const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
832  RCP<t_vector_t> tv = xtv->getTpetra_Vector();
833 
834  RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
835  *tv, numLocalRows, myNewRows);
836 
837  RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
838 
839  return xvnew;
840  }
841  }
842 };
843 
845 // Tpetra::MultiVector
846 template <typename scalar_t,
847  typename lno_t,
848  typename gno_t,
849  typename node_t>
850 struct XpetraTraits<Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
851 {
852  typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
853  typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
854  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
855 
856  static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
857  {
858  return rcp(new xt_vector_t(a));
859  }
860 
861  static RCP<t_vector_t> doMigration(const t_vector_t &from,
862  size_t numLocalElts, const gno_t *myNewElts)
863  {
864  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
865 
866  // source map
867  const RCP<const map_t> &smap = from.getMap();
868  gno_t numGlobalElts = smap->getGlobalNumElements();
869  gno_t base = smap->getMinAllGlobalIndex();
870 
871  // target map
872  ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
873  const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
874  RCP<const map_t> tmap = rcp(new map_t(numGlobalElts, eltList, base, comm));
875 
876  // importer
877  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
878 
879  // target vector
880  RCP<t_vector_t> MV = rcp(
881  new t_vector_t(tmap, from.getNumVectors(), true));
882  MV->doImport(from, importer, Tpetra::INSERT);
883 
884  return MV;
885  }
886 };
887 
889 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
890 // Epetra_MultiVector
891 template < >
892 struct XpetraTraits<Epetra_MultiVector>
893 {
894  typedef InputTraits<Epetra_MultiVector>::lno_t lno_t;
895  typedef InputTraits<Epetra_MultiVector>::gno_t gno_t;
896  typedef InputTraits<Epetra_MultiVector>::node_t node_t;
897  typedef InputTraits<Epetra_MultiVector>::scalar_t scalar_t;
898  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
899 
900  static inline RCP<x_mvector_t> convertToXpetra(
901  const RCP<Epetra_MultiVector> &a)
902  {
903  RCP<Xpetra::EpetraMultiVectorT<gno_t, node_t> > xemv;
904  try {
905  xemv = rcp(new Xpetra::EpetraMultiVectorT<gno_t,node_t>(a));
906  }
907  catch (std::exception &e) {
908  if (std::is_same<node_t, Xpetra::EpetraNode>::value)
909  throw std::runtime_error(std::string("Cannot convert from "
910  "Epetra_MultiVector to "
911  "Xpetra::EpetraMultiVectorT\n")
912  + e.what());
913  else
914  throw std::runtime_error(std::string("Cannot convert from "
915  "Epetra_MultiVector to "
916  "Xpetra::EpetraMultiVectorT\n"
917  "Use node_t that is supported by "
918  "Xpetra with Epetra classes\n")
919  + e.what());
920  }
921  return rcp_implicit_cast<x_mvector_t>(xemv);
922  }
923 
924  static RCP<Epetra_MultiVector> doMigration(const Epetra_MultiVector &from,
925  size_t numLocalElts, const gno_t *myNewElts)
926  {
927  // source map
928  const Epetra_BlockMap &smap = from.Map();
929  gno_t numGlobalElts = smap.NumGlobalElements();
930  int base = smap.MinAllGID();
931 
932  // target map
933  const Epetra_Comm &comm = from.Comm();
934  const Epetra_BlockMap tmap(numGlobalElts, numLocalElts, myNewElts,
935  1, base, comm);
936 
937  // importer
938  Epetra_Import importer(tmap, smap);
939 
940  // target vector
941  RCP<Epetra_MultiVector> MV = rcp(
942  new Epetra_MultiVector(tmap, from.NumVectors(), true));
943  Epetra_CombineMode c = Insert;
944  MV->Import(from, importer, c);
945 
946  return MV;
947  }
948 };
949 #endif
950 
952 // Xpetra::MultiVector
953 template <typename scalar_t,
954  typename lno_t,
955  typename gno_t,
956  typename node_t>
957 struct XpetraTraits<Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
958 {
959  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
960  typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
961  typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
962 
963  static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
964  {
965  return a;
966  }
967 
968  static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
969  size_t numLocalRows, const gno_t *myNewRows)
970  {
971  Xpetra::UnderlyingLib lib = from.getMap()->lib();
972 
973  if (lib == Xpetra::UseEpetra){
974  throw std::logic_error("compiler should have used specialization");
975  } else{
976  // Do the import with the Tpetra::MultiVector traits object
977  const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
978  RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
979 
980  RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
981  *tv, numLocalRows, myNewRows);
982 
983  RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
984 
985  return xvnew;
986  }
987  }
988 };
989 
991 // Xpetra::MultiVector specialization
992 template <typename node_t>
993 struct XpetraTraits<Xpetra::MultiVector<double, int, int, node_t> >
994 {
995  typedef double scalar_t;
996  typedef int lno_t;
997  typedef int gno_t;
998  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
999  typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
1000  typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
1001 
1002  static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
1003  {
1004  return a;
1005  }
1006 
1007  static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
1008  size_t numLocalRows, const gno_t *myNewRows)
1009  {
1010  Xpetra::UnderlyingLib lib = from.getMap()->lib();
1011 
1012  if (lib == Xpetra::UseEpetra){
1013 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
1014  typedef Epetra_MultiVector e_mvector_t;
1015  typedef Xpetra::EpetraMultiVectorT<gno_t,node_t> xe_mvector_t;
1016  // Do the import with the Epetra_MultiVector traits object
1017  const xe_mvector_t *xev = dynamic_cast<const xe_mvector_t *>(&from);
1018  RCP<e_mvector_t> ev = xev->getEpetra_MultiVector();
1019 
1020  RCP<e_mvector_t> evnew = XpetraTraits<e_mvector_t>::doMigration(
1021  *ev, numLocalRows, myNewRows);
1022 
1023  RCP<x_mvector_t> xvnew = XpetraTraits<e_mvector_t>::convertToXpetra(evnew);
1024 
1025  return xvnew;
1026 #else
1027  throw std::runtime_error("Xpetra with Epetra requested, but "
1028  "Trilinos is not built with Epetra");
1029 #endif
1030  } else{
1031  // Do the import with the Tpetra::MultiVector traits object
1032  const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
1033  RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
1034 
1035  RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
1036  *tv, numLocalRows, myNewRows);
1037 
1038  RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
1039 
1040  return xvnew;
1041  }
1042  }
1043 };
1044 
1045 #endif // DOXYGEN_SHOULD_SKIP_THIS
1046 
1047 } //namespace Zoltan2
1048 
1049 #endif // _ZOLTAN2_XPETRATRAITS_HPP_
Defines the traits required for Tpetra, Eptra and Xpetra objects.
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
default_gno_t gno_t
The objects global ordinal data type.
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tmatrix_t
Xpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > xmatrix_t
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...
Xpetra::CrsGraph< zlno_t, zgno_t, znode_t > xgraph_t
::Tpetra::Details::DefaultTypes::global_ordinal_type default_gno_t
::Tpetra::Details::DefaultTypes::local_ordinal_type default_lno_t
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
dictionary vals
Definition: xml2dox.py:186
Traits for application input objects.
default_lno_t lno_t
The objects local ordinal data type.
Gathering definitions used in software development.
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tgraph_t