Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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_DefaultComm.hpp>
58 #include <Teuchos_RCP.hpp>
59 #include <Teuchos_Array.hpp>
60 #include <Teuchos_ArrayRCP.hpp>
61 #include <Teuchos_Comm.hpp>
62 #include <Teuchos_VerboseObject.hpp>
63 #include <Tpetra_CrsMatrix.hpp>
64 #include <Tpetra_Vector.hpp>
65 
66 #ifdef HAVE_ZOLTAN_EPETRA
67 #include <Xpetra_EpetraUtils.hpp>
68 #endif
69 
70 using std::string;
71 using Teuchos::RCP;
72 using Teuchos::ArrayRCP;
73 using Teuchos::ArrayView;
74 using Teuchos::Array;
75 using Teuchos::rcp;
76 using Teuchos::Comm;
77 
78 ArrayRCP<zgno_t> roundRobinMapShared(
79  int proc,
80  int nprocs,
81  zgno_t basegid,
82  zgno_t maxgid,
83  size_t nglobalrows
84 )
85 {
86  if (nglobalrows != size_t(maxgid - basegid + 1)){
87  std::cout << "Error: Map is invalid for test - fix test" << std::endl;
88  std::cerr << "Error: Map is invalid for test - fix test" << std::endl;
89  std::cout << "FAIL" << std::endl;
90  exit(1);
91  }
92  RCP<Array<zgno_t> > mygids = rcp(new Array<zgno_t>);
93  zgno_t firstzgno_t = proc;
94  if (firstzgno_t < basegid){
95  zgno_t n = basegid % proc;
96  if (n>0)
97  firstzgno_t = basegid - n + proc;
98  else
99  firstzgno_t = basegid;
100  }
101  for (zgno_t gid=firstzgno_t; gid <= maxgid; gid+=nprocs){
102  (*mygids).append(gid);
103  }
104 
105  ArrayRCP<zgno_t> newIdArcp = Teuchos::arcp(mygids);
106 
107  return newIdArcp;
108 }
109 
110 #ifdef HAVE_EPETRA_DATA_TYPES
111 ArrayRCP<zgno_t> roundRobinMap(const Epetra_BlockMap &emap)
112 {
113  const Epetra_Comm &comm = emap.Comm();
114  int proc = comm.MyPID();
115  int nprocs = comm.NumProc();
116  zgno_t basegid = emap.MinAllGID();
117  zgno_t maxgid = emap.MaxAllGID();
118  size_t nglobalrows = emap.NumGlobalElements();
119 
120  return roundRobinMapShared(proc, nprocs, basegid, maxgid, nglobalrows);
121 }
122 #endif
123 
124 ArrayRCP<zgno_t> roundRobinMap(const Tpetra::Map<zlno_t, zgno_t, znode_t> &tmap)
125 {
126  const RCP<const Comm<int> > &comm = tmap.getComm();
127  int proc = comm->getRank();
128  int nprocs = comm->getSize();
129  zgno_t basegid = tmap.getMinAllGlobalIndex();
130  zgno_t maxgid = tmap.getMaxAllGlobalIndex();
131  size_t nglobalrows = tmap.getGlobalNumElements();
132 
133  return roundRobinMapShared(proc, nprocs, basegid, maxgid, nglobalrows);
134 }
135 
136 ArrayRCP<zgno_t> roundRobinMap(const Xpetra::Map<zlno_t, zgno_t, znode_t> &xmap)
137 {
138  const RCP<const Comm<int> > &comm = xmap.getComm();
139  int proc = comm->getRank();
140  int nprocs = comm->getSize();
141  zgno_t basegid = xmap.getMinAllGlobalIndex();
142  zgno_t maxgid = xmap.getMaxAllGlobalIndex();
143  size_t nglobalrows = xmap.getGlobalNumElements();
144 
145  return roundRobinMapShared(proc, nprocs, basegid, maxgid, nglobalrows);
146 }
147 
148 int main(int narg, char *arg[])
149 {
150  Tpetra::ScopeGuard tscope(&narg, &arg);
151  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
152 
153  int rank = comm->getRank();
154  bool aok = true;
155 
156  Teuchos::RCP<Teuchos::FancyOStream> outStream =
157  Teuchos::VerboseObjectBase::getDefaultOStream();
158  Teuchos::EVerbosityLevel v=Teuchos::VERB_EXTREME;
159 
160  typedef Tpetra::CrsMatrix<zscalar_t,zlno_t,zgno_t,znode_t> tmatrix_t;
161  typedef Tpetra::CrsGraph<zlno_t,zgno_t,znode_t> tgraph_t;
162  typedef Tpetra::Vector<zscalar_t,zlno_t,zgno_t,znode_t> tvector_t;
163  typedef Tpetra::MultiVector<zscalar_t,zlno_t,zgno_t,znode_t> tmvector_t;
164  typedef Xpetra::CrsMatrix<zscalar_t,zlno_t,zgno_t,znode_t> xmatrix_t;
165  typedef Xpetra::CrsGraph<zlno_t,zgno_t,znode_t> xgraph_t;
166  typedef Xpetra::Vector<zscalar_t,zlno_t,zgno_t,znode_t> xvector_t;
167  typedef Xpetra::MultiVector<zscalar_t,zlno_t,zgno_t,znode_t> xmvector_t;
168 
169  // Create object that can give us test Tpetra and Xpetra input.
170 
171  RCP<UserInputForTests> uinput;
172  Teuchos::ParameterList params;
173  params.set("input file", "simple");
174  params.set("file type", "Chaco");
175 
176  try{
177  uinput = rcp(new UserInputForTests(params, comm));
178  }
179  catch(std::exception &e){
180  aok = false;
181  std::cout << e.what() << std::endl;
182  }
183  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
184 
186  // Tpetra::CrsMatrix
187  // Tpetra::CrsGraph
188  // Tpetra::Vector
189  // Tpetra::MultiVector
191 
192 
193  // XpetraTraits<Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> >
194  {
195  RCP<tmatrix_t> M;
196 
197  try{
198  M = uinput->getUITpetraCrsMatrix();
199  }
200  catch(std::exception &e){
201  aok = false;
202  std::cout << e.what() << std::endl;
203  }
204  TEST_FAIL_AND_EXIT(*comm, aok, "getTpetraCrsMatrix ", 1);
205 
206  if (rank== 0)
207  std::cout << "Original Tpetra matrix " << M->getGlobalNumRows()
208  << " x " << M->getGlobalNumCols() << std::endl;
209 
210  M->describe(*outStream,v);
211 
212  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(M->getRowMap()));
213 
214  zgno_t localNumRows = newRowIds.size();
215 
216  RCP<const tmatrix_t> newM;
217  try{
219  localNumRows, newRowIds.getRawPtr());
220  }
221  catch(std::exception &e){
222  aok = false;
223  std::cout << e.what() << std::endl;
224  }
225  TEST_FAIL_AND_EXIT(*comm, aok,
226  " Zoltan2::XpetraTraits<tmatrix_t>::doMigration ", 1);
227 
228  if (rank== 0)
229  std::cout << "Migrated Tpetra matrix" << std::endl;
230 
231  newM->describe(*outStream,v);
232  }
233 
234  // XpetraTraits<Tpetra::CrsGraph<zscalar_t, zlno_t, zgno_t, znode_t> >
235  {
236  RCP<tgraph_t> G;
237 
238  try{
239  G = uinput->getUITpetraCrsGraph();
240  }
241  catch(std::exception &e){
242  aok = false;
243  std::cout << e.what() << std::endl;
244  }
245  TEST_FAIL_AND_EXIT(*comm, aok, "getTpetraCrsGraph ", 1);
246 
247  if (rank== 0)
248  std::cout << "Original Tpetra graph" << std::endl;
249 
250  G->describe(*outStream,v);
251 
252  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(G->getRowMap()));
253 
254  zgno_t localNumRows = newRowIds.size();
255 
256  RCP<const tgraph_t> newG;
257  try{
259  localNumRows, newRowIds.getRawPtr());
260  }
261  catch(std::exception &e){
262  aok = false;
263  std::cout << e.what() << std::endl;
264  }
265  TEST_FAIL_AND_EXIT(*comm, aok,
266  " Zoltan2::XpetraTraits<tgraph_t>::doMigration ", 1);
267 
268  if (rank== 0)
269  std::cout << "Migrated Tpetra graph" << std::endl;
270 
271  newG->describe(*outStream,v);
272  }
273 
274  // XpetraTraits<Tpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>>
275  {
276  RCP<tvector_t> V;
277 
278  try{
279  V = rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 1));
280  V->randomize();
281  }
282  catch(std::exception &e){
283  aok = false;
284  std::cout << e.what() << std::endl;
285  }
286  TEST_FAIL_AND_EXIT(*comm, aok, "getTpetraVector", 1);
287 
288  if (rank== 0)
289  std::cout << "Original Tpetra vector" << std::endl;
290 
291  V->describe(*outStream,v);
292 
293  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(V->getMap()));
294 
295  zgno_t localNumRows = newRowIds.size();
296 
297  RCP<const tvector_t> newV;
298  try{
300  localNumRows, newRowIds.getRawPtr());
301  }
302  catch(std::exception &e){
303  aok = false;
304  std::cout << e.what() << std::endl;
305  }
306  TEST_FAIL_AND_EXIT(*comm, aok,
307  " Zoltan2::XpetraTraits<tvector_t>::doMigration ", 1);
308 
309  if (rank== 0)
310  std::cout << "Migrated Tpetra vector" << std::endl;
311 
312  newV->describe(*outStream,v);
313  }
314 
315  // XpetraTraits<Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>>
316  {
317  RCP<tmvector_t> MV;
318 
319  try{
320  MV = rcp(new tmvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 3));
321  MV->randomize();
322  }
323  catch(std::exception &e){
324  aok = false;
325  std::cout << e.what() << std::endl;
326  }
327  TEST_FAIL_AND_EXIT(*comm, aok, "getTpetraMultiVector", 1);
328 
329  if (rank== 0)
330  std::cout << "Original Tpetra multivector" << std::endl;
331 
332  MV->describe(*outStream,v);
333 
334  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(MV->getMap()));
335 
336  zgno_t localNumRows = newRowIds.size();
337 
338  RCP<const tmvector_t> newMV;
339  try{
341  localNumRows, newRowIds.getRawPtr());
342  }
343  catch(std::exception &e){
344  aok = false;
345  std::cout << e.what() << std::endl;
346  }
347  TEST_FAIL_AND_EXIT(*comm, aok,
348  " Zoltan2::XpetraTraits<tmvector_t>::doMigration ", 1);
349 
350  if (rank== 0)
351  std::cout << "Migrated Tpetra multivector" << std::endl;
352 
353  newMV->describe(*outStream,v);
354  }
355 
357  // Xpetra::CrsMatrix
358  // Xpetra::CrsGraph
359  // Xpetra::Vector
360  // Xpetra::MultiVector
362 
363  // XpetraTraits<Xpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> >
364  {
365  RCP<xmatrix_t> M;
366 
367  try{
368  M = uinput->getUIXpetraCrsMatrix();
369  }
370  catch(std::exception &e){
371  aok = false;
372  std::cout << e.what() << std::endl;
373  }
374  TEST_FAIL_AND_EXIT(*comm, aok, "getXpetraCrsMatrix ", 1);
375 
376  if (rank== 0)
377  std::cout << "Original Xpetra matrix" << std::endl;
378 
379  M->describe(*outStream,v);
380 
381  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(M->getRowMap()));
382 
383  zgno_t localNumRows = newRowIds.size();
384 
385  RCP<const xmatrix_t> newM;
386  try{
388  localNumRows, newRowIds.getRawPtr());
389  }
390  catch(std::exception &e){
391  aok = false;
392  std::cout << e.what() << std::endl;
393  }
394  TEST_FAIL_AND_EXIT(*comm, aok,
395  " Zoltan2::XpetraTraits<xmatrix_t>::doMigration ", 1);
396 
397  if (rank== 0)
398  std::cout << "Migrated Xpetra matrix" << std::endl;
399 
400  newM->describe(*outStream,v);
401  }
402 
403  // XpetraTraits<Xpetra::CrsGraph<zscalar_t, zlno_t, zgno_t, znode_t> >
404  {
405  RCP<xgraph_t> G;
406 
407  try{
408  G = uinput->getUIXpetraCrsGraph();
409  }
410  catch(std::exception &e){
411  aok = false;
412  std::cout << e.what() << std::endl;
413  }
414  TEST_FAIL_AND_EXIT(*comm, aok, "getXpetraCrsGraph ", 1);
415 
416  if (rank== 0)
417  std::cout << "Original Xpetra graph" << std::endl;
418 
419  G->describe(*outStream,v);
420 
421  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(G->getRowMap()));
422 
423  zgno_t localNumRows = newRowIds.size();
424 
425  RCP<const xgraph_t> newG;
426  try{
428  localNumRows, newRowIds.getRawPtr());
429  }
430  catch(std::exception &e){
431  aok = false;
432  std::cout << e.what() << std::endl;
433  }
434  TEST_FAIL_AND_EXIT(*comm, aok,
435  " Zoltan2::XpetraTraits<xgraph_t>::doMigration ", 1);
436 
437  if (rank== 0)
438  std::cout << "Migrated Xpetra graph" << std::endl;
439 
440  newG->describe(*outStream,v);
441  }
442 
443  // XpetraTraits<Xpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>>
444  {
445  RCP<xvector_t> V;
446 
447  try{
448  RCP<tvector_t> tV =
449  rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 1));
450  tV->randomize();
452  }
453  catch(std::exception &e){
454  aok = false;
455  std::cout << e.what() << std::endl;
456  }
457  TEST_FAIL_AND_EXIT(*comm, aok, "getXpetraVector", 1);
458 
459  if (rank== 0)
460  std::cout << "Original Xpetra vector" << std::endl;
461 
462  V->describe(*outStream,v);
463 
464  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(V->getMap()));
465 
466  zgno_t localNumRows = newRowIds.size();
467 
468  RCP<const xvector_t> newV;
469  try{
471  localNumRows, newRowIds.getRawPtr());
472  }
473  catch(std::exception &e){
474  aok = false;
475  std::cout << e.what() << std::endl;
476  }
477  TEST_FAIL_AND_EXIT(*comm, aok,
478  " Zoltan2::XpetraTraits<xvector_t>::doMigration ", 1);
479 
480  if (rank== 0)
481  std::cout << "Migrated Xpetra vector" << std::endl;
482 
483  newV->describe(*outStream,v);
484  }
485 
486  // XpetraTraits<Xpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>>
487  {
488  RCP<xmvector_t> MV;
489 
490  try{
491  RCP<tmvector_t> tMV =
492  rcp(new tmvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 3));
493  tMV->randomize();
495  }
496  catch(std::exception &e){
497  aok = false;
498  std::cout << e.what() << std::endl;
499  }
500  TEST_FAIL_AND_EXIT(*comm, aok, "getXpetraMultiVector", 1);
501 
502  if (rank== 0)
503  std::cout << "Original Xpetra multivector" << std::endl;
504 
505  MV->describe(*outStream,v);
506 
507  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(MV->getMap()));
508 
509  zgno_t localNumRows = newRowIds.size();
510 
511  RCP<const xmvector_t> newMV;
512  try{
514  localNumRows, newRowIds.getRawPtr());
515  }
516  catch(std::exception &e){
517  aok = false;
518  std::cout << e.what() << std::endl;
519  }
520  TEST_FAIL_AND_EXIT(*comm, aok,
521  " Zoltan2::XpetraTraits<xmvector_t>::doMigration ", 1);
522 
523  if (rank== 0)
524  std::cout << "Migrated Xpetra multivector" << std::endl;
525 
526  newMV->describe(*outStream,v);
527  }
528 
529 #ifdef HAVE_EPETRA_DATA_TYPES
530  // Epetra_CrsMatrix
532  // Epetra_CrsGraph
533  // Epetra_Vector
534  // Epetra_MultiVector
536 
537  typedef Epetra_CrsMatrix ematrix_t;
538  typedef Epetra_CrsGraph egraph_t;
539  typedef Epetra_Vector evector_t;
540  typedef Epetra_MultiVector emvector_t;
541  typedef Epetra_BlockMap emap_t;
542 
543  // XpetraTraits<Epetra_CrsMatrix>
544  {
545  RCP<ematrix_t> M;
546 
547  try{
548  M = uinput->getUIEpetraCrsMatrix();
549  }
550  catch(std::exception &e){
551  aok = false;
552  std::cout << e.what() << std::endl;
553  }
554  TEST_FAIL_AND_EXIT(*comm, aok, "getEpetraCrsMatrix ", 1);
555 
556  if (rank== 0)
557  std::cout << "Original Epetra matrix" << std::endl;
558 
559  M->Print(std::cout);
560 
561  RCP<const emap_t> emap = Teuchos::rcpFromRef(M->RowMap());
562  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*emap);
563 
564  zgno_t localNumRows = newRowIds.size();
565 
566  RCP<const ematrix_t> newM;
567  try{
569  localNumRows, newRowIds.getRawPtr());
570  }
571  catch(std::exception &e){
572  aok = false;
573  std::cout << e.what() << std::endl;
574  }
575  TEST_FAIL_AND_EXIT(*comm, aok,
576  " Zoltan2::XpetraTraits<ematrix_t>::doMigration ", 1);
577 
578  if (rank== 0)
579  std::cout << "Migrated Epetra matrix" << std::endl;
580 
581  newM->Print(std::cout);
582  }
583 
584  // XpetraTraits<Epetra_CrsGraph>
585  {
586  RCP<egraph_t> G;
587 
588  try{
589  G = uinput->getUIEpetraCrsGraph();
590  }
591  catch(std::exception &e){
592  aok = false;
593  std::cout << e.what() << std::endl;
594  }
595  TEST_FAIL_AND_EXIT(*comm, aok, "getEpetraCrsGraph ", 1);
596 
597  if (rank== 0)
598  std::cout << "Original Epetra graph" << std::endl;
599 
600  G->Print(std::cout);
601 
602  RCP<const emap_t> emap = Teuchos::rcpFromRef(G->RowMap());
603  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*emap);
604 
605  zgno_t localNumRows = newRowIds.size();
606 
607  RCP<const egraph_t> newG;
608  try{
610  localNumRows, newRowIds.getRawPtr());
611  }
612  catch(std::exception &e){
613  aok = false;
614  std::cout << e.what() << std::endl;
615  }
616  TEST_FAIL_AND_EXIT(*comm, aok,
617  " Zoltan2::XpetraTraits<egraph_t>::doMigration ", 1);
618 
619  if (rank== 0)
620  std::cout << "Migrated Epetra graph" << std::endl;
621 
622  newG->Print(std::cout);
623  }
624 
625  // XpetraTraits<Epetra_Vector>
626  {
627  RCP<evector_t> V;
628 
629  try{
630  V = rcp(new Epetra_Vector(uinput->getUIEpetraCrsGraph()->RowMap()));
631  V->Random();
632  }
633  catch(std::exception &e){
634  aok = false;
635  std::cout << e.what() << std::endl;
636  }
637  TEST_FAIL_AND_EXIT(*comm, aok, "getEpetraVector", 1);
638 
639  if (rank== 0)
640  std::cout << "Original Epetra vector" << std::endl;
641 
642  V->Print(std::cout);
643 
644  RCP<const emap_t> emap = Teuchos::rcpFromRef(V->Map());
645  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*emap);
646 
647  zgno_t localNumRows = newRowIds.size();
648 
649  RCP<const evector_t> newV;
650  try{
652  localNumRows, newRowIds.getRawPtr());
653  }
654  catch(std::exception &e){
655  aok = false;
656  std::cout << e.what() << std::endl;
657  }
658  TEST_FAIL_AND_EXIT(*comm, aok,
659  " Zoltan2::XpetraTraits<evector_t>::doMigration ", 1);
660 
661  if (rank== 0)
662  std::cout << "Migrated Epetra vector" << std::endl;
663 
664  newV->Print(std::cout);
665  }
666 
667  // XpetraTraits<Epetra_MultiVector>
668  {
669  RCP<emvector_t> MV;
670 
671  try{
672  MV =
673  rcp(new Epetra_MultiVector(uinput->getUIEpetraCrsGraph()->RowMap(),3));
674  MV->Random();
675  }
676  catch(std::exception &e){
677  aok = false;
678  std::cout << e.what() << std::endl;
679  }
680  TEST_FAIL_AND_EXIT(*comm, aok, "getEpetraMultiVector", 1);
681 
682  if (rank== 0)
683  std::cout << "Original Epetra multivector" << std::endl;
684 
685  MV->Print(std::cout);
686 
687  RCP<const emap_t> emap = Teuchos::rcpFromRef(MV->Map());
688  ArrayRCP<zgno_t> newRowIds = roundRobinMap(*emap);
689 
690  zgno_t localNumRows = newRowIds.size();
691 
692  RCP<const emvector_t> newMV;
693  try{
695  localNumRows, newRowIds.getRawPtr());
696  }
697  catch(std::exception &e){
698  aok = false;
699  std::cout << e.what() << std::endl;
700  }
701  TEST_FAIL_AND_EXIT(*comm, aok,
702  " Zoltan2::XpetraTraits<emvector_t>::doMigration ", 1);
703 
704  if (rank== 0)
705  std::cout << "Migrated Epetra multivector" << std::endl;
706 
707  newMV->Print(std::cout);
708  }
709 #endif // have epetra data types (int, int, double)
710 
712  // DONE
714 
715  if (rank==0)
716  std::cout << "PASS" << std::endl;
717 }
718 
ArrayRCP< zgno_t > roundRobinMapShared(int proc, int nprocs, zgno_t basegid, zgno_t maxgid, size_t nglobalrows)
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.
int main(int narg, char **arg)
Definition: coloring1.cpp:199
common code used by tests
ArrayRCP< zgno_t > roundRobinMap(const Tpetra::Map< zlno_t, zgno_t, znode_t > &tmap)
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
Xpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
Tpetra::CrsGraph< zlno_t, zgno_t, znode_t > tgraph_t
Tpetra::Map::global_ordinal_type zgno_t