Zoltan2
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_EpetraCrsMatrix.hpp>
57 #include <Xpetra_TpetraCrsMatrix.hpp>
58 #include <Xpetra_TpetraRowMatrix.hpp>
59 #include <Xpetra_EpetraVector.hpp>
60 #include <Xpetra_TpetraVector.hpp>
61 #include <Xpetra_EpetraUtils.hpp>
62 #include <Tpetra_Vector.hpp>
63 
64 namespace Zoltan2 {
65 
67 // Extra traits needed only for Xpetra matrices and graphs
68 
91 template <typename User>
93 {
96  static inline RCP<User> convertToXpetra(const RCP<User> &a)
97  {
98  return a;
99  }
100 
104 
108 
114  static RCP<User> doMigration(const User &from,
115  size_t numLocalRows, const gno_t *myNewRows)
116  {
117  return Teuchos::null;
118  }
119 };
120 
121 #ifndef DOXYGEN_SHOULD_SKIP_THIS
122 
124 // Tpetra::CrsMatrix
125 template <typename scalar_t,
126  typename lno_t,
127  typename gno_t,
128  typename node_t>
129 struct XpetraTraits<Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
130 {
131  typedef typename Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t>
132  xmatrix_t;
133  typedef typename Xpetra::TpetraCrsMatrix<scalar_t,lno_t,gno_t,node_t>
134  xtmatrix_t;
135  typedef typename Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t>
136  tmatrix_t;
137 
138  static inline RCP<xmatrix_t> convertToXpetra(const RCP<tmatrix_t> &a)
139  {
140  return rcp(new xtmatrix_t(a));
141  }
142 
143  static RCP<tmatrix_t> doMigration(const tmatrix_t &from,
144  size_t numLocalRows, const gno_t *myNewRows)
145  {
146  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
147  lno_t base = 0;
148 
149  // source map
150  const RCP<const map_t> &smap = from.getRowMap();
151  int oldNumElts = smap->getNodeNumElements();
152  gno_t numGlobalRows = smap->getGlobalNumElements();
153 
154  // target map
155  ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
156  const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
157  RCP<const map_t> tmap = rcp(
158  new map_t(numGlobalRows, rowList, base, comm));
159  int newNumElts = numLocalRows;
160 
161  // importer
162  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
163 
164  // number of non zeros in my new rows
165  typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
166  vector_t numOld(smap); // TODO These vectors should have scalar = size_t,
167  vector_t numNew(tmap); // but explicit instantiation does not yet support that.
168  for (int lid=0; lid < oldNumElts; lid++){
169  numOld.replaceGlobalValue(smap->getGlobalElement(lid),
170  scalar_t(from.getNumEntriesInLocalRow(lid)));
171  }
172  numNew.doImport(numOld, importer, Tpetra::INSERT);
173 
174  // TODO Could skip this copy if could declare vector with scalar = size_t.
175  ArrayRCP<size_t> nnz(newNumElts);
176  if (newNumElts > 0){
177  ArrayRCP<scalar_t> ptr = numNew.getDataNonConst(0);
178  for (int lid=0; lid < newNumElts; lid++){
179  nnz[lid] = static_cast<size_t>(ptr[lid]);
180  }
181  }
182 
183  // target matrix
184  RCP<tmatrix_t> M = rcp(new tmatrix_t(tmap, nnz, Tpetra::StaticProfile));
185  M->doImport(from, importer, Tpetra::INSERT);
186  M->fillComplete();
187 
188  return M;
189  }
190 };
191 
193 // Epetra_CrsMatrix
194 
195 template <>
196 struct XpetraTraits<Epetra_CrsMatrix>
197 {
202 
203  static inline RCP<Xpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> >
204  convertToXpetra(const RCP<Epetra_CrsMatrix> &a)
205  {
206  return rcp(new Xpetra::EpetraCrsMatrix(a));
207  }
208 
209 
210  static RCP<Epetra_CrsMatrix> doMigration(const Epetra_CrsMatrix &from,
211  size_t numLocalRows, const gno_t *myNewRows)
212  {
213  lno_t base = 0;
214 
215  // source map
216  const Epetra_Map &smap = from.RowMap();
217  int oldNumElts = smap.NumMyElements();
218  gno_t numGlobalRows = smap.NumGlobalElements();
219 
220  // target map
221  const Epetra_Comm &comm = from.Comm();
222  Epetra_Map tmap(numGlobalRows, numLocalRows, myNewRows, base, comm);
223  int newNumElts = numLocalRows;
224 
225  // importer
226  Epetra_Import importer(tmap, smap);
227 
228  // number of non zeros in my new rows
229  Epetra_Vector numOld(smap);
230  Epetra_Vector numNew(tmap);
231 
232  for (int lid=0; lid < oldNumElts; lid++){
233  numOld[lid] = from.NumMyEntries(lid);
234  }
235  numNew.Import(numOld, importer, Insert);
236 
237  Array<int> nnz(newNumElts);
238  for (int lid=0; lid < newNumElts; lid++){
239  nnz[lid] = static_cast<int>(numNew[lid]);
240  }
241 
242  // target matrix
243  RCP<Epetra_CrsMatrix> M = rcp(new Epetra_CrsMatrix(::Copy, tmap, nnz.getRawPtr(), true));
244  M->Import(from, importer, Insert);
245  M->FillComplete();
246 
247  return M;
248  }
249 };
250 
252 // Xpetra::CrsMatrix
253 // KDDKDD: Do we need specializations for Xpetra::EpetraCrsMatrix and
254 // KDDKDD: Xpetra::TpetraCrsMatrix
255 template <typename scalar_t,
256  typename lno_t,
257  typename gno_t,
258  typename node_t>
259 struct XpetraTraits<Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> >
260 {
261  typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
262  typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
263  typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
264 
265  static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
266  {
267  return a;
268  }
269 
270  static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
271  size_t numLocalRows, const gno_t *myNewRows)
272  {
273  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
274 
275  if (lib == Xpetra::UseEpetra){
276  throw std::logic_error("compiler should have used specialization");
277  } else{
278  // Do the import with the Tpetra::CrsMatrix traits object
279  const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
280  RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
281 
282  RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
283  *tm, numLocalRows, myNewRows);
284 
285  RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
286 
287  return xmnew;
288  }
289  }
290 };
291 
293 // Xpetra::CrsMatrix specialization
294 
295 template <typename node_t>
296 struct XpetraTraits<Xpetra::CrsMatrix<double, int, int, node_t> >
297 {
298  typedef double scalar_t;
299  typedef int lno_t;
300  typedef int gno_t;
301  typedef Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
302  typedef Xpetra::TpetraCrsMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
303  typedef Tpetra::CrsMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
304  typedef Xpetra::EpetraCrsMatrix xe_matrix_t;
305  typedef Epetra_CrsMatrix e_matrix_t;
306 
307  static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
308  {
309  return a;
310  }
311 
312  static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
313  size_t numLocalRows, const gno_t *myNewRows)
314  {
315  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
316 
317  if (lib == Xpetra::UseEpetra){
318  // Do the import with the Epetra_CrsMatrix traits object
319  const xe_matrix_t *xem = dynamic_cast<const xe_matrix_t *>(&from);
320  RCP<const e_matrix_t> em = xem->getEpetra_CrsMatrix();
321 
322  RCP<e_matrix_t> emnew = XpetraTraits<e_matrix_t>::doMigration(
323  *em, numLocalRows, myNewRows);
324 
325  RCP<x_matrix_t> xmnew = XpetraTraits<e_matrix_t>::convertToXpetra(emnew);
326 
327  return xmnew;
328  } else{
329  // Do the import with the Tpetra::CrsMatrix traits object
330  const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
331  RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
332 
333  RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
334  *tm, numLocalRows, myNewRows);
335 
336  RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
337 
338  return xmnew;
339  }
340  }
341 };
342 
343 
345 // Tpetra::CrsGraph
346 template <typename lno_t,
347  typename gno_t,
348  typename node_t>
349 struct XpetraTraits<Tpetra::CrsGraph<lno_t, gno_t, node_t> >
350 {
351  typedef typename Xpetra::CrsGraph<lno_t, gno_t, node_t> xgraph_t;
352  typedef typename Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xtgraph_t;
353  typedef typename Tpetra::CrsGraph<lno_t, gno_t, node_t> tgraph_t;
354 
355  static inline RCP<xgraph_t> convertToXpetra(const RCP<tgraph_t> &a)
356  {
357  return rcp(new xtgraph_t(a));
358  }
359 
360  static RCP<tgraph_t> doMigration(const tgraph_t &from,
361  size_t numLocalRows, const gno_t *myNewRows)
362  {
363  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
364  lno_t base = 0;
365 
366  // source map
367  const RCP<const map_t> &smap = from.getRowMap();
368  int oldNumElts = smap->getNodeNumElements();
369  gno_t numGlobalRows = smap->getGlobalNumElements();
370 
371  // target map
372  ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
373  const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
374  RCP<const map_t> tmap = rcp(
375  new map_t(numGlobalRows, rowList, base, comm));
376 
377  // importer
378  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
379 
380  // number of entries in my new rows
381  typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
382  vector_t numOld(smap);
383  vector_t numNew(tmap);
384  for (int lid=0; lid < oldNumElts; lid++){
385  numOld.replaceGlobalValue(smap->getGlobalElement(lid),
386  from.getNumEntriesInLocalRow(lid));
387  }
388  numNew.doImport(numOld, importer, Tpetra::INSERT);
389 
390  size_t numElts = tmap->getNodeNumElements();
391  ArrayRCP<const gno_t> nnz;
392  if (numElts > 0)
393  nnz = numNew.getData(0); // hangs if vector len == 0
394 
395  ArrayRCP<const size_t> nnz_size_t;
396 
397  if (numElts && sizeof(gno_t) != sizeof(size_t)){
398  size_t *vals = new size_t [numElts];
399  nnz_size_t = arcp(vals, 0, numElts, true);
400  for (size_t i=0; i < numElts; i++){
401  vals[i] = static_cast<size_t>(nnz[i]);
402  }
403  }
404  else{
405  nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
406  }
407 
408  // target graph
409  RCP<tgraph_t> G = rcp(new tgraph_t(tmap, nnz_size_t, Tpetra::StaticProfile));
410 
411  G->doImport(from, importer, Tpetra::INSERT);
412  G->fillComplete();
413 
414  return G;
415  }
416 
417 };
418 
420 // Epetra_CrsGraph
421 template < >
422 struct XpetraTraits<Epetra_CrsGraph>
423 {
427  static inline RCP<Xpetra::CrsGraph<lno_t,gno_t,node_t> >
428  convertToXpetra(const RCP<Epetra_CrsGraph> &a)
429  {
430  return rcp(new Xpetra::EpetraCrsGraph(a));
431  }
432 
433  static RCP<Epetra_CrsGraph> doMigration(const Epetra_CrsGraph &from,
434  size_t numLocalRows, const gno_t *myNewRows)
435  {
436  lno_t base = 0;
437 
438  // source map
439  const Epetra_BlockMap &smap = from.RowMap();
440  gno_t numGlobalRows = smap.NumGlobalElements();
441  lno_t oldNumElts = smap.NumMyElements();
442 
443  // target map
444  const Epetra_Comm &comm = from.Comm();
445  Epetra_BlockMap tmap(numGlobalRows, numLocalRows,
446  myNewRows, 1, base, comm);
447  lno_t newNumElts = tmap.NumMyElements();
448 
449  // importer
450  Epetra_Import importer(tmap, smap);
451 
452  // number of non zeros in my new rows
453  Epetra_Vector numOld(smap);
454  Epetra_Vector numNew(tmap);
455 
456  for (int lid=0; lid < oldNumElts; lid++){
457  numOld[lid] = from.NumMyIndices(lid);
458  }
459  numNew.Import(numOld, importer, Insert);
460 
461  Array<int> nnz(newNumElts);
462  for (int lid=0; lid < newNumElts; lid++){
463  nnz[lid] = static_cast<int>(numNew[lid]);
464  }
465 
466  // target graph
467  RCP<Epetra_CrsGraph> G = rcp(new Epetra_CrsGraph(::Copy, tmap, nnz.getRawPtr(), true));
468  G->Import(from, importer, Insert);
469  G->FillComplete();
470 
471  return G;
472  }
473 
474 };
475 
477 // Xpetra::CrsGraph
478 // KDDKDD Do we need specializations for Xpetra::TpetraCrsGraph and
479 // KDDKDD Xpetra::EpetraCrsGraph?
480 template <typename lno_t,
481  typename gno_t,
482  typename node_t>
483 struct XpetraTraits<Xpetra::CrsGraph<lno_t, gno_t, node_t> >
484 {
485  typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
486  typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
487  typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
488 
489  static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
490  {
491  return a;
492  }
493 
494  static RCP<x_graph_t> doMigration(const x_graph_t &from,
495  size_t numLocalRows, const gno_t *myNewRows)
496  {
497  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
498 
499  if (lib == Xpetra::UseEpetra){
500  throw std::logic_error("compiler should have used specialization");
501  } else{
502  // Do the import with the Tpetra::CrsGraph traits object
503  const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
504  RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
505 
506  RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
507  *tg, numLocalRows, myNewRows);
508 
509  RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
510  return xgnew;
511  }
512  }
513 };
514 
515 
516 
518 // Xpetra::RowMatrix
519 template <typename scalar_t,
520  typename lno_t,
521  typename gno_t,
522  typename node_t>
523 struct XpetraTraits<Xpetra::RowMatrix<scalar_t, lno_t, gno_t, node_t> >
524 {
525  typedef Xpetra::RowMatrix<scalar_t, lno_t, gno_t, node_t> x_matrix_t;
526  typedef Xpetra::TpetraRowMatrix<scalar_t, lno_t, gno_t, node_t> xt_matrix_t;
527  typedef Tpetra::RowMatrix<scalar_t,lno_t,gno_t,node_t> t_matrix_t;
528 
529  static inline RCP<x_matrix_t> convertToXpetra(const RCP<x_matrix_t > &a)
530  {
531  return a;
532  }
533 
534  static RCP<x_matrix_t> doMigration(const x_matrix_t &from,
535  size_t numLocalRows, const gno_t *myNewRows)
536  {
537  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
538 
539  if (lib == Xpetra::UseEpetra){
540  throw std::logic_error("compiler should have used specialization");
541  } else{
542  // Do the import with the Tpetra::CrsMatrix traits object
543  const xt_matrix_t *xtm = dynamic_cast<const xt_matrix_t *>(&from);
544  RCP<const t_matrix_t> tm = xtm->getTpetra_CrsMatrix();
545 
546  RCP<t_matrix_t> tmnew = XpetraTraits<t_matrix_t>::doMigration(
547  *tm, numLocalRows, myNewRows);
548 
549  RCP<x_matrix_t> xmnew = XpetraTraits<t_matrix_t>::convertToXpetra(tmnew);
550 
551  return xmnew;
552  }
553  }
554 };
555 
557 // Xpetra::CrsGraph specialization
558 template < typename node_t>
559 struct XpetraTraits<Xpetra::CrsGraph<int, int, node_t> >
560 {
561  typedef int lno_t;
562  typedef int gno_t;
563  typedef Xpetra::CrsGraph<lno_t, gno_t, node_t> x_graph_t;
564  typedef Xpetra::TpetraCrsGraph<lno_t, gno_t, node_t> xt_graph_t;
565  typedef Tpetra::CrsGraph<lno_t,gno_t,node_t> t_graph_t;
566  typedef Xpetra::EpetraCrsGraph xe_graph_t;
567  typedef Epetra_CrsGraph e_graph_t;
568 
569  static inline RCP<x_graph_t> convertToXpetra(const RCP<x_graph_t> &a)
570  {
571  return a;
572  }
573 
574  static RCP<x_graph_t> doMigration(const x_graph_t &from,
575  size_t numLocalRows, const gno_t *myNewRows)
576  {
577  Xpetra::UnderlyingLib lib = from.getRowMap()->lib();
578 
579  if (lib == Xpetra::UseEpetra){
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  // Do the import with the Tpetra::CrsGraph traits object
592  const xt_graph_t *xtg = dynamic_cast<const xt_graph_t *>(&from);
593  RCP<const t_graph_t> tg = xtg->getTpetra_CrsGraph();
594 
595  RCP<t_graph_t> tgnew = XpetraTraits<t_graph_t>::doMigration(
596  *tg, numLocalRows, myNewRows);
597 
598  RCP<x_graph_t> xgnew = XpetraTraits<t_graph_t>::convertToXpetra(tgnew);
599 
600  return xgnew;
601  }
602  }
603 };
604 
606 // Tpetra::Vector
607 template <typename scalar_t,
608  typename lno_t,
609  typename gno_t,
610  typename node_t>
611 struct XpetraTraits<Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
612 {
613  typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
614  typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
615  typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
616 
617  static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
618  {
619  return rcp(new xt_vector_t(a));
620  }
621 
622  static RCP<t_vector_t> doMigration(const t_vector_t &from,
623  size_t numLocalElts, const gno_t *myNewElts)
624  {
625  lno_t base = 0;
626  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
627 
628  // source map
629  const RCP<const map_t> &smap = from.getMap();
630  gno_t numGlobalElts = smap->getGlobalNumElements();
631 
632  // target map
633  ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
634  const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
635  RCP<const map_t> tmap = rcp(
636  new map_t(numGlobalElts, eltList, base, comm));
637 
638  // importer
639  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
640 
641  // target vector
642  RCP<t_vector_t> V =
643  Tpetra::createVector<scalar_t,lno_t,gno_t,node_t>(tmap);
644  V->doImport(from, importer, Tpetra::INSERT);
645 
646  return V;
647  }
648 };
649 
651 // Epetra_Vector
652 template < >
653 struct XpetraTraits<Epetra_Vector>
654 {
657  typedef InputTraits<Epetra_Vector>::node_t node_t;
658  typedef InputTraits<Epetra_Vector>::scalar_t scalar_t;
659 
660  typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
661 
662  static inline RCP<x_vector_t> convertToXpetra(const RCP<Epetra_Vector> &a)
663  {
664  RCP<Xpetra::EpetraVector> xev = rcp(new Xpetra::EpetraVector(a));
665  return rcp_implicit_cast<x_vector_t>(xev);
666  }
667 
668  static RCP<Epetra_Vector> doMigration(const Epetra_Vector &from,
669  size_t numLocalElts, const gno_t *myNewElts)
670  {
671  lno_t base = 0;
672  // source map
673  const Epetra_BlockMap &smap = from.Map();
674  gno_t numGlobalElts = smap.NumGlobalElements();
675 
676  // target map
677  const Epetra_Comm &comm = from.Comm();
678  const Epetra_BlockMap tmap(numGlobalElts, numLocalElts, myNewElts,
679  1, base, comm);
680 
681  // importer
682  Epetra_Import importer(tmap, smap);
683 
684  // target vector
685  RCP<Epetra_Vector> V = rcp(new Epetra_Vector(tmap, true));
686  Epetra_CombineMode c = Insert;
687  V->Import(from, importer, c);
688 
689  return V;
690  }
691 };
692 
694 // Xpetra::Vector
695 template <typename scalar_t,
696  typename lno_t,
697  typename gno_t,
698  typename node_t>
699 struct XpetraTraits<Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> >
700 {
701  typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
702  typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
703  typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
704 
705  static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
706  {
707  return a;
708  }
709 
710  static RCP<x_vector_t> doMigration(const x_vector_t &from,
711  size_t numLocalRows, const gno_t *myNewRows)
712  {
713  Xpetra::UnderlyingLib lib = from.getMap()->lib();
714 
715  if (lib == Xpetra::UseEpetra){
716  throw std::logic_error("compiler should have used specialization");
717  } else{
718  // Do the import with the Tpetra::Vector traits object
719  const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
720  RCP<const t_vector_t> tv = xtv->getTpetra_Vector();
721 
722  RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
723  *tv, numLocalRows, myNewRows);
724 
725  RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
726 
727  return xvnew;
728  }
729  }
730 };
731 
733 // Xpetra::Vector specialization
734 template <typename node_t>
735 struct XpetraTraits<Xpetra::Vector<double, int, int, node_t> >
736 {
737  typedef double scalar_t;
738  typedef int lno_t;
739  typedef int gno_t;
740  typedef Xpetra::Vector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
741  typedef Xpetra::TpetraVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
742  typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
743  typedef Xpetra::EpetraVector xe_vector_t;
744  typedef Epetra_Vector e_vector_t;
745 
746  static inline RCP<x_vector_t> convertToXpetra(const RCP<x_vector_t> &a)
747  {
748  return a;
749  }
750 
751  static RCP<x_vector_t> doMigration(const x_vector_t &from,
752  size_t numLocalRows, const gno_t *myNewRows)
753  {
754  Xpetra::UnderlyingLib lib = from.getMap()->lib();
755 
756  if (lib == Xpetra::UseEpetra){
757  // Do the import with the Epetra_Vector traits object
758  const xe_vector_t *xev = dynamic_cast<const xe_vector_t *>(&from);
759  RCP<const e_vector_t> ev = rcp(xev->getEpetra_Vector());
760 
761  RCP<e_vector_t> evnew = XpetraTraits<e_vector_t>::doMigration(
762  *ev, numLocalRows, myNewRows);
763 
764  RCP<x_vector_t> xvnew = XpetraTraits<e_vector_t>::convertToXpetra(evnew);
765 
766  return xvnew;
767  } else{
768  // Do the import with the Tpetra::Vector traits object
769  const xt_vector_t *xtv = dynamic_cast<const xt_vector_t *>(&from);
770  RCP<t_vector_t> tv = xtv->getTpetra_Vector();
771 
772  RCP<t_vector_t> tvnew = XpetraTraits<t_vector_t>::doMigration(
773  *tv, numLocalRows, myNewRows);
774 
775  RCP<x_vector_t> xvnew = XpetraTraits<t_vector_t>::convertToXpetra(tvnew);
776 
777  return xvnew;
778  }
779  }
780 };
781 
783 // Tpetra::MultiVector
784 template <typename scalar_t,
785  typename lno_t,
786  typename gno_t,
787  typename node_t>
788 struct XpetraTraits<Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
789 {
790  typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_vector_t;
791  typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_vector_t;
792  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_vector_t;
793 
794  static inline RCP<x_vector_t> convertToXpetra(const RCP<t_vector_t> &a)
795  {
796  return rcp(new xt_vector_t(a));
797  }
798 
799  static RCP<t_vector_t> doMigration(const t_vector_t &from,
800  size_t numLocalElts, const gno_t *myNewElts)
801  {
802  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
803  lno_t base = 0;
804 
805  // source map
806  const RCP<const map_t> &smap = from.getMap();
807  gno_t numGlobalElts = smap->getGlobalNumElements();
808 
809  // target map
810  ArrayView<const gno_t> eltList(myNewElts, numLocalElts);
811  const RCP<const Teuchos::Comm<int> > comm = from.getMap()->getComm();
812  RCP<const map_t> tmap = rcp(
813  new map_t(numGlobalElts, eltList, base, comm));
814 
815  // importer
816  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
817 
818  // target vector
819  RCP<t_vector_t> MV = rcp(
820  new t_vector_t(tmap, from.getNumVectors(), true));
821  MV->doImport(from, importer, Tpetra::INSERT);
822 
823  return MV;
824  }
825 };
826 
828 // Epetra_MultiVector
829 template < >
830 struct XpetraTraits<Epetra_MultiVector>
831 {
836  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
837 
838  static inline RCP<x_mvector_t> convertToXpetra(
839  const RCP<Epetra_MultiVector> &a)
840  {
841  RCP<Xpetra::EpetraMultiVector> xemv = rcp(new Xpetra::EpetraMultiVector(a));
842  return rcp_implicit_cast<x_mvector_t>(xemv);
843  }
844 
845  static RCP<Epetra_MultiVector> doMigration(const Epetra_MultiVector &from,
846  size_t numLocalElts, const gno_t *myNewElts)
847  {
848  lno_t base = 0;
849  // source map
850  const Epetra_BlockMap &smap = from.Map();
851  gno_t numGlobalElts = smap.NumGlobalElements();
852 
853  // target map
854  const Epetra_Comm &comm = from.Comm();
855  const Epetra_BlockMap tmap(numGlobalElts, numLocalElts, myNewElts,
856  1, base, comm);
857 
858  // importer
859  Epetra_Import importer(tmap, smap);
860 
861  // target vector
862  RCP<Epetra_MultiVector> MV = rcp(
863  new Epetra_MultiVector(tmap, from.NumVectors(), true));
864  Epetra_CombineMode c = Insert;
865  MV->Import(from, importer, c);
866 
867  return MV;
868  }
869 };
870 
872 // Xpetra::MultiVector
873 template <typename scalar_t,
874  typename lno_t,
875  typename gno_t,
876  typename node_t>
877 struct XpetraTraits<Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> >
878 {
879  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
880  typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
881  typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
882 
883  static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
884  {
885  return a;
886  }
887 
888  static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
889  size_t numLocalRows, const gno_t *myNewRows)
890  {
891  Xpetra::UnderlyingLib lib = from.getMap()->lib();
892 
893  if (lib == Xpetra::UseEpetra){
894  throw std::logic_error("compiler should have used specialization");
895  } else{
896  // Do the import with the Tpetra::MultiVector traits object
897  const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
898  RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
899 
900  RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
901  *tv, numLocalRows, myNewRows);
902 
903  RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
904 
905  return xvnew;
906  }
907  }
908 };
909 
911 // Xpetra::MultiVector specialization
912 template <typename node_t>
913 struct XpetraTraits<Xpetra::MultiVector<double, int, int, node_t> >
914 {
915  typedef double scalar_t;
916  typedef int lno_t;
917  typedef int gno_t;
918  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
919  typedef Xpetra::TpetraMultiVector<scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
920  typedef Tpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> t_mvector_t;
921  typedef Xpetra::EpetraMultiVector xe_mvector_t;
922  typedef Epetra_MultiVector e_mvector_t;
923 
924  static inline RCP<x_mvector_t> convertToXpetra(const RCP<x_mvector_t> &a)
925  {
926  return a;
927  }
928 
929  static RCP<x_mvector_t> doMigration(const x_mvector_t &from,
930  size_t numLocalRows, const gno_t *myNewRows)
931  {
932  Xpetra::UnderlyingLib lib = from.getMap()->lib();
933 
934  if (lib == Xpetra::UseEpetra){
935  // Do the import with the Epetra_MultiVector traits object
936  const xe_mvector_t *xev = dynamic_cast<const xe_mvector_t *>(&from);
937  RCP<e_mvector_t> ev = xev->getEpetra_MultiVector();
938 
939  RCP<e_mvector_t> evnew = XpetraTraits<e_mvector_t>::doMigration(
940  *ev, numLocalRows, myNewRows);
941 
942  RCP<x_mvector_t> xvnew = XpetraTraits<e_mvector_t>::convertToXpetra(evnew);
943 
944  return xvnew;
945 
946  } else{
947  // Do the import with the Tpetra::MultiVector traits object
948  const xt_mvector_t *xtv = dynamic_cast<const xt_mvector_t *>(&from);
949  RCP<t_mvector_t> tv = xtv->getTpetra_MultiVector();
950 
951  RCP<t_mvector_t> tvnew = XpetraTraits<t_mvector_t>::doMigration(
952  *tv, numLocalRows, myNewRows);
953 
954  RCP<x_mvector_t> xvnew = XpetraTraits<t_mvector_t>::convertToXpetra(tvnew);
955 
956  return xvnew;
957  }
958  }
959 };
960 
961 #endif // DOXYGEN_SHOULD_SKIP_THIS
962 
963 } //namespace Zoltan2
964 
965 #endif // _ZOLTAN2_XPETRATRAITS_HPP_
Defines the traits required for Tpetra, Eptra and Xpetra objects.
default_gno_t gno_t
The objects global ordinal data type.
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
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
dictionary vals
Definition: xml2dox.py:186
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices...
Traits for application input objects.
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 objects local ordinal data type.
Tpetra::Vector< z2TestScalar, z2TestLO, z2TestGO > Vector
Gathering definitions used in software development.
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tgraph_t
default_scalar_t scalar_t
The data type for weights and coordinates.