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