Zoltan2
XpetraTraits.cpp
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 //
46 // Basic test of the XpetraTraits definitions.
47 //
48 // TODO - a real test would figure out if the migrated objects are
49 // the same as the original, here we just look at them on stdout.
50 // TODO look at number of diagonals and max number of entries in
51 // Tpetra and Xpetra migrated graphs. They're garbage.
52 
53 #include <Zoltan2_XpetraTraits.hpp>
54 #include <Zoltan2_TestHelpers.hpp>
55 
56 #include <string>
57 #include <Teuchos_GlobalMPISession.hpp>
58 #include <Teuchos_DefaultComm.hpp>
59 #include <Teuchos_RCP.hpp>
60 #include <Teuchos_Array.hpp>
61 #include <Teuchos_ArrayRCP.hpp>
62 #include <Teuchos_Comm.hpp>
63 #include <Teuchos_VerboseObject.hpp>
64 #include <Tpetra_CrsMatrix.hpp>
65 #include <Tpetra_Vector.hpp>
66 
67 #include <Xpetra_EpetraUtils.hpp>
68 #ifdef HAVE_ZOLTAN2_MPI
69 #include <Epetra_MpiComm.h>
70 #else
71 #include <Epetra_SerialComm.h>
72 #endif
73 
74 using namespace std;
75 using std::string;
76 using Teuchos::RCP;
77 using Teuchos::ArrayRCP;
78 using Teuchos::ArrayView;
79 using Teuchos::Array;
80 using Teuchos::rcp;
81 using Teuchos::Comm;
82 
83 ArrayRCP<zgno_t> roundRobinMap(
84  const RCP<const Xpetra::Map<zlno_t, zgno_t, znode_t> > &m)
85 {
86  const RCP<const Comm<int> > &comm = m->getComm();
87  int proc = comm->getRank();
88  int nprocs = comm->getSize();
89  zgno_t base = m->getMinAllGlobalIndex();
90  zgno_t max = m->getMaxAllGlobalIndex();
91  size_t globalrows = m->getGlobalNumElements();
92  if (globalrows != size_t(max - base + 1)){
93  TEST_FAIL_AND_EXIT(*comm, 0,
94  string("Map is invalid for test - fix test"), 1);
95  }
96  RCP<Array<zgno_t> > mygids = rcp(new Array<zgno_t>);
97  zgno_t firstzgno_t = proc;
98  if (firstzgno_t < base){
99  zgno_t n = base % proc;
100  if (n>0)
101  firstzgno_t = base - n + proc;
102  else
103  firstzgno_t = base;
104  }
105  for (zgno_t gid=firstzgno_t; gid <= max; gid+=nprocs){
106  (*mygids).append(gid);
107  }
108 
109  ArrayRCP<zgno_t> newIdArcp = Teuchos::arcp(mygids);
110 
111  return newIdArcp;
112 }
113 
114 int main(int argc, char *argv[])
115 {
116  Teuchos::GlobalMPISession session(&argc, &argv);
117  RCP<const Comm<int> > comm = Teuchos::DefaultComm<int>::getComm();
118  int rank = comm->getRank();
119 
120  Teuchos::RCP<Teuchos::FancyOStream> outStream =
121  Teuchos::VerboseObjectBase::getDefaultOStream();
122  Teuchos::EVerbosityLevel v=Teuchos::VERB_EXTREME;
123 
124  typedef Tpetra::CrsMatrix<zscalar_t,zlno_t,zgno_t,znode_t> tmatrix_t;
125  typedef Tpetra::CrsGraph<zlno_t,zgno_t,znode_t> tgraph_t;
126  typedef Tpetra::Vector<zscalar_t,zlno_t,zgno_t,znode_t> tvector_t;
127  typedef Tpetra::MultiVector<zscalar_t,zlno_t,zgno_t,znode_t> tmvector_t;
128  typedef Xpetra::CrsMatrix<zscalar_t,zlno_t,zgno_t,znode_t> xmatrix_t;
129  typedef Xpetra::CrsGraph<zlno_t,zgno_t,znode_t> xgraph_t;
130  typedef Xpetra::Vector<zscalar_t,zlno_t,zgno_t,znode_t> xvector_t;
131  typedef Xpetra::MultiVector<zscalar_t,zlno_t,zgno_t,znode_t> xmvector_t;
132  typedef Xpetra::TpetraMap<zlno_t,zgno_t,znode_t> xtmap_t;
133 
134  // Create object that can give us test Tpetra and Xpetra input.
135 
136  RCP<UserInputForTests> uinput;
137 
138  try{
139  uinput =
140  rcp(new UserInputForTests(testDataFilePath,std::string("simple"), comm, true));
141  }
142  catch(std::exception &e){
143  TEST_FAIL_AND_EXIT(*comm, 0, string("input ")+e.what(), 1);
144  }
145 
147  // Tpetra::CrsMatrix
148  // Tpetra::CrsGraph
149  // Tpetra::Vector
150  // Tpetra::MultiVector
152 
153  // XpetraTraits<Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> >
154  {
155  RCP<tmatrix_t> M;
156 
157  try{
158  M = uinput->getUITpetraCrsMatrix();
159  }
160  catch(std::exception &e){
161  TEST_FAIL_AND_EXIT(*comm, 0,
162  string("getTpetraCrsMatrix ")+e.what(), 1);
163  }
164 
165  if (rank== 0)
166  std::cout << "Original Tpetra matrix " << M->getGlobalNumRows()
167  << " x " << M->getGlobalNumCols() << std::endl;
168 
169  M->describe(*outStream,v);
170 
171  RCP<const xtmap_t> xmap(new xtmap_t(M->getRowMap()));
172 
173  ArrayRCP<zgno_t> newRowIds = roundRobinMap(xmap);
174 
175  zgno_t localNumRows = newRowIds.size();
176 
177  RCP<const tmatrix_t> newM;
178  try{
180  localNumRows, newRowIds.getRawPtr());
181  }
182  catch(std::exception &e){
183  TEST_FAIL_AND_EXIT(*comm, 0,
184  string(" Zoltan2::XpetraTraits<tmatrix_t>::doMigration ")+e.what(), 1);
185  }
186 
187  if (rank== 0)
188  std::cout << "Migrated Tpetra matrix" << std::endl;
189 
190  newM->describe(*outStream,v);
191  }
192 
193  // XpetraTraits<Tpetra::CrsGraph<zscalar_t, zlno_t, zgno_t, znode_t> >
194  {
195  RCP<tgraph_t> G;
196 
197  try{
198  G = uinput->getUITpetraCrsGraph();
199  }
200  catch(std::exception &e){
201  TEST_FAIL_AND_EXIT(*comm, 0,
202  string("getTpetraCrsGraph ")+e.what(), 1);
203  }
204 
205  if (rank== 0)
206  std::cout << "Original Tpetra graph" << std::endl;
207 
208  G->describe(*outStream,v);
209 
210  RCP<const xtmap_t> xmap(new xtmap_t(G->getRowMap()));
211  ArrayRCP<zgno_t> newRowIds = roundRobinMap(xmap);
212 
213  zgno_t localNumRows = newRowIds.size();
214 
215  RCP<const tgraph_t> newG;
216  try{
218  localNumRows, newRowIds.getRawPtr());
219  }
220  catch(std::exception &e){
221  TEST_FAIL_AND_EXIT(*comm, 0,
222  string(" Zoltan2::XpetraTraits<tgraph_t>::doMigration ")+e.what(), 1);
223  }
224 
225  if (rank== 0)
226  std::cout << "Migrated Tpetra graph" << std::endl;
227 
228  newG->describe(*outStream,v);
229  }
230 
231  // XpetraTraits<Tpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>>
232  {
233  RCP<tvector_t> V;
234 
235  try{
236  V = rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 1));
237  V->randomize();
238  }
239  catch(std::exception &e){
240  TEST_FAIL_AND_EXIT(*comm, 0,
241  string("getTpetraVector")+e.what(), 1);
242  }
243 
244  if (rank== 0)
245  std::cout << "Original Tpetra vector" << std::endl;
246 
247  V->describe(*outStream,v);
248 
249  RCP<const xtmap_t> xmap(new xtmap_t(V->getMap()));
250  ArrayRCP<zgno_t> newRowIds = roundRobinMap(xmap);
251 
252  zgno_t localNumRows = newRowIds.size();
253 
254  RCP<const tvector_t> newV;
255  try{
257  localNumRows, newRowIds.getRawPtr());
258  }
259  catch(std::exception &e){
260  TEST_FAIL_AND_EXIT(*comm, 0,
261  string(" Zoltan2::XpetraTraits<tvector_t>::doMigration ")+e.what(), 1);
262  }
263 
264  if (rank== 0)
265  std::cout << "Migrated Tpetra vector" << std::endl;
266 
267  newV->describe(*outStream,v);
268  }
269 
270  // XpetraTraits<Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>>
271  {
272  RCP<tmvector_t> MV;
273 
274  try{
275  MV = rcp(new tmvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 3));
276  MV->randomize();
277  }
278  catch(std::exception &e){
279  TEST_FAIL_AND_EXIT(*comm, 0,
280  string("getTpetraMultiVector")+e.what(), 1);
281  }
282 
283  if (rank== 0)
284  std::cout << "Original Tpetra multivector" << std::endl;
285 
286  MV->describe(*outStream,v);
287 
288  RCP<const xtmap_t> xmap(new xtmap_t(MV->getMap()));
289  ArrayRCP<zgno_t> newRowIds = roundRobinMap(xmap);
290 
291  zgno_t localNumRows = newRowIds.size();
292 
293  RCP<const tmvector_t> newMV;
294  try{
296  localNumRows, newRowIds.getRawPtr());
297  }
298  catch(std::exception &e){
299  TEST_FAIL_AND_EXIT(*comm, 0,
300  string(" Zoltan2::XpetraTraits<tmvector_t>::doMigration ")+e.what(), 1);
301  }
302 
303  if (rank== 0)
304  std::cout << "Migrated Tpetra multivector" << std::endl;
305 
306  newMV->describe(*outStream,v);
307  }
308 
310  // Xpetra::CrsMatrix
311  // Xpetra::CrsGraph
312  // Xpetra::Vector
313  // Xpetra::MultiVector
315 
316  // XpetraTraits<Xpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> >
317  {
318  RCP<xmatrix_t> M;
319 
320  try{
321  M = uinput->getUIXpetraCrsMatrix();
322  }
323  catch(std::exception &e){
324  TEST_FAIL_AND_EXIT(*comm, 0,
325  string("getXpetraCrsMatrix ")+e.what(), 1);
326  }
327 
328  if (rank== 0)
329  std::cout << "Original Xpetra matrix" << std::endl;
330 
331  M->describe(*outStream,v);
332 
333  ArrayRCP<zgno_t> newRowIds = roundRobinMap(M->getRowMap());
334 
335  zgno_t localNumRows = newRowIds.size();
336 
337  RCP<const xmatrix_t> newM;
338  try{
340  localNumRows, newRowIds.getRawPtr());
341  }
342  catch(std::exception &e){
343  TEST_FAIL_AND_EXIT(*comm, 0,
344  string(" Zoltan2::XpetraTraits<xmatrix_t>::doMigration ")+e.what(), 1);
345  }
346 
347  if (rank== 0)
348  std::cout << "Migrated Xpetra matrix" << std::endl;
349 
350  newM->describe(*outStream,v);
351  }
352 
353  // XpetraTraits<Xpetra::CrsGraph<zscalar_t, zlno_t, zgno_t, znode_t> >
354  {
355  RCP<xgraph_t> G;
356 
357  try{
358  G = uinput->getUIXpetraCrsGraph();
359  }
360  catch(std::exception &e){
361  TEST_FAIL_AND_EXIT(*comm, 0,
362  string("getXpetraCrsGraph ")+e.what(), 1);
363  }
364 
365  if (rank== 0)
366  std::cout << "Original Xpetra graph" << std::endl;
367 
368  G->describe(*outStream,v);
369 
370  ArrayRCP<zgno_t> newRowIds = roundRobinMap(G->getRowMap());
371 
372  zgno_t localNumRows = newRowIds.size();
373 
374  RCP<const xgraph_t> newG;
375  try{
377  localNumRows, newRowIds.getRawPtr());
378  }
379  catch(std::exception &e){
380  TEST_FAIL_AND_EXIT(*comm, 0,
381  string(" Zoltan2::XpetraTraits<xgraph_t>::doMigration ")+e.what(), 1);
382  }
383 
384  if (rank== 0)
385  std::cout << "Migrated Xpetra graph" << std::endl;
386 
387  newG->describe(*outStream,v);
388  }
389 
390  // XpetraTraits<Xpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>>
391  {
392  RCP<xvector_t> V;
393 
394  try{
395  RCP<tvector_t> tV =
396  rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 1));
397  tV->randomize();
399  }
400  catch(std::exception &e){
401  TEST_FAIL_AND_EXIT(*comm, 0,
402  string("getXpetraVector")+e.what(), 1);
403  }
404 
405  if (rank== 0)
406  std::cout << "Original Xpetra vector" << std::endl;
407 
408  V->describe(*outStream,v);
409 
410  ArrayRCP<zgno_t> newRowIds = roundRobinMap(V->getMap());
411 
412  zgno_t localNumRows = newRowIds.size();
413 
414  RCP<const xvector_t> newV;
415  try{
417  localNumRows, newRowIds.getRawPtr());
418  }
419  catch(std::exception &e){
420  TEST_FAIL_AND_EXIT(*comm, 0,
421  string(" Zoltan2::XpetraTraits<xvector_t>::doMigration ")+e.what(), 1);
422  }
423 
424  if (rank== 0)
425  std::cout << "Migrated Xpetra vector" << std::endl;
426 
427  newV->describe(*outStream,v);
428  }
429 
430  // XpetraTraits<Xpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>>
431  {
432  RCP<xmvector_t> MV;
433 
434  try{
435  RCP<tmvector_t> tMV =
436  rcp(new tmvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 3));
437  tMV->randomize();
439  }
440  catch(std::exception &e){
441  TEST_FAIL_AND_EXIT(*comm, 0,
442  string("getXpetraMultiVector")+e.what(), 1);
443  }
444 
445  if (rank== 0)
446  std::cout << "Original Xpetra multivector" << std::endl;
447 
448  MV->describe(*outStream,v);
449 
450  ArrayRCP<zgno_t> newRowIds = roundRobinMap(MV->getMap());
451 
452  zgno_t localNumRows = newRowIds.size();
453 
454  RCP<const xmvector_t> newMV;
455  try{
457  localNumRows, newRowIds.getRawPtr());
458  }
459  catch(std::exception &e){
460  TEST_FAIL_AND_EXIT(*comm, 0,
461  string(" Zoltan2::XpetraTraits<xmvector_t>::doMigration ")+e.what(), 1);
462  }
463 
464  if (rank== 0)
465  std::cout << "Migrated Xpetra multivector" << std::endl;
466 
467  newMV->describe(*outStream,v);
468  }
469 
470 #ifdef HAVE_EPETRA_DATA_TYPES
471  // Epetra_CrsMatrix
473  // Epetra_CrsGraph
474  // Epetra_Vector
475  // Epetra_MultiVector
477 
478  typedef Epetra_CrsMatrix ematrix_t;
479  typedef Epetra_CrsGraph egraph_t;
480  typedef Epetra_Vector evector_t;
481  typedef Epetra_MultiVector emvector_t;
482  typedef Xpetra::EpetraMap xemap_t;
483  typedef Epetra_BlockMap emap_t;
484 
485  // Create object that can give us test Epetra input.
486 
487  RCP<UserInputForTests> euinput;
488 
489  try{
490  euinput =
491  rcp(new UserInputForTests(testDataFilePath,std::string("simple"), comm, true));
492  }
493  catch(std::exception &e){
494  TEST_FAIL_AND_EXIT(*comm, 0, string("epetra input ")+e.what(), 1);
495  }
496 
497  // XpetraTraits<Epetra_CrsMatrix>
498  {
499  RCP<ematrix_t> M;
500 
501  try{
502  M = euinput->getUIEpetraCrsMatrix();
503  }
504  catch(std::exception &e){
505  TEST_FAIL_AND_EXIT(*comm, 0,
506  string("getEpetraCrsMatrix ")+e.what(), 1);
507  }
508 
509  if (rank== 0)
510  std::cout << "Original Epetra matrix" << std::endl;
511 
512  M->Print(std::cout);
513 
514  RCP<const emap_t> emap = Teuchos::rcpFromRef(M->RowMap());
515  RCP<const xemap_t> xmap(new xemap_t(emap));
516 
517  ArrayRCP<zgno_t> newRowIds = roundRobinMap(xmap);
518 
519  zgno_t localNumRows = newRowIds.size();
520 
521  RCP<const ematrix_t> newM;
522  try{
524  localNumRows, newRowIds.getRawPtr());
525  }
526  catch(std::exception &e){
527  TEST_FAIL_AND_EXIT(*comm, 0,
528  string(" Zoltan2::XpetraTraits<ematrix_t>::doMigration ")+e.what(), 1);
529  }
530 
531  if (rank== 0)
532  std::cout << "Migrated Epetra matrix" << std::endl;
533 
534  newM->Print(std::cout);
535  }
536 
537  // XpetraTraits<Epetra_CrsGraph>
538  {
539  RCP<egraph_t> G;
540 
541  try{
542  G = euinput->getUIEpetraCrsGraph();
543  }
544  catch(std::exception &e){
545  TEST_FAIL_AND_EXIT(*comm, 0,
546  string("getEpetraCrsGraph ")+e.what(), 1);
547  }
548 
549  if (rank== 0)
550  std::cout << "Original Epetra graph" << std::endl;
551 
552  G->Print(std::cout);
553 
554  RCP<const emap_t> emap = Teuchos::rcpFromRef(G->RowMap());
555  RCP<const xemap_t> xmap(new xemap_t(emap));
556  ArrayRCP<zgno_t> newRowIds = roundRobinMap(xmap);
557 
558  zgno_t localNumRows = newRowIds.size();
559 
560  RCP<const egraph_t> newG;
561  try{
563  localNumRows, newRowIds.getRawPtr());
564  }
565  catch(std::exception &e){
566  TEST_FAIL_AND_EXIT(*comm, 0,
567  string(" Zoltan2::XpetraTraits<egraph_t>::doMigration ")+e.what(), 1);
568  }
569 
570  if (rank== 0)
571  std::cout << "Migrated Epetra graph" << std::endl;
572 
573  newG->Print(std::cout);
574  }
575 
576  // XpetraTraits<Epetra_Vector>
577  {
578  RCP<evector_t> V;
579 
580  try{
581  V = rcp(new Epetra_Vector(euinput->getUIEpetraCrsGraph()->RowMap()));
582  V->Random();
583  }
584  catch(std::exception &e){
585  TEST_FAIL_AND_EXIT(*comm, 0,
586  string("getEpetraVector")+e.what(), 1);
587  }
588 
589  if (rank== 0)
590  std::cout << "Original Epetra vector" << std::endl;
591 
592  V->Print(std::cout);
593 
594  RCP<const emap_t> emap = Teuchos::rcpFromRef(V->Map());
595  RCP<const xemap_t> xmap(new xemap_t(emap));
596  ArrayRCP<zgno_t> newRowIds = roundRobinMap(xmap);
597 
598  zgno_t localNumRows = newRowIds.size();
599 
600  RCP<const evector_t> newV;
601  try{
603  localNumRows, newRowIds.getRawPtr());
604  }
605  catch(std::exception &e){
606  TEST_FAIL_AND_EXIT(*comm, 0,
607  string(" Zoltan2::XpetraTraits<evector_t>::doMigration ")+e.what(), 1);
608  }
609 
610  if (rank== 0)
611  std::cout << "Migrated Epetra vector" << std::endl;
612 
613  newV->Print(std::cout);
614  }
615 
616  // XpetraTraits<Epetra_MultiVector>
617  {
618  RCP<emvector_t> MV;
619 
620  try{
621  MV =
622  rcp(new Epetra_MultiVector(euinput->getUIEpetraCrsGraph()->RowMap(),3));
623  MV->Random();
624  }
625  catch(std::exception &e){
626  TEST_FAIL_AND_EXIT(*comm, 0,
627  string("getEpetraMultiVector")+e.what(), 1);
628  }
629 
630  if (rank== 0)
631  std::cout << "Original Epetra multivector" << std::endl;
632 
633  MV->Print(std::cout);
634 
635  RCP<const emap_t> emap = Teuchos::rcpFromRef(MV->Map());
636  RCP<const xemap_t> xmap(new xemap_t(emap));
637  ArrayRCP<zgno_t> newRowIds = roundRobinMap(xmap);
638 
639  zgno_t localNumRows = newRowIds.size();
640 
641  RCP<const emvector_t> newMV;
642  try{
644  localNumRows, newRowIds.getRawPtr());
645  }
646  catch(std::exception &e){
647  TEST_FAIL_AND_EXIT(*comm, 0,
648  string(" Zoltan2::XpetraTraits<emvector_t>::doMigration ")+e.what(), 1);
649  }
650 
651  if (rank== 0)
652  std::cout << "Migrated Epetra multivector" << std::endl;
653 
654  newMV->Print(std::cout);
655  }
656 #endif // have epetra data types (int, int, double)
657 
659  // DONE
661 
662  if (rank==0)
663  std::cout << "PASS" << std::endl;
664 }
665 
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tmatrix_t
Xpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > xmatrix_t
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
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...
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
Xpetra::CrsGraph< zlno_t, zgno_t, znode_t > xgraph_t
Traits of Xpetra classes, including migration method.
common code used by tests
Epetra_CrsGraph egraph_t
int main(int argc, char *argv[])
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
ArrayRCP< zgno_t > roundRobinMap(const RCP< const Xpetra::Map< zlno_t, zgno_t, znode_t > > &m)
Epetra_MultiVector evector_t
int zgno_t
Xpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
Epetra_CrsMatrix ematrix_t
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tgraph_t
std::string testDataFilePath(".")