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
31using std::string;
32using Teuchos::RCP;
33using Teuchos::ArrayRCP;
34using Teuchos::ArrayView;
35using Teuchos::Array;
36using Teuchos::rcp;
37using Teuchos::Comm;
38
39ArrayRCP<zgno_t> roundRobinMapShared(
40 int proc,
41 int nprocs,
42 zgno_t basegid,
43 zgno_t maxgid,
44 size_t nglobalrows
45)
46{
47 if (nglobalrows != size_t(maxgid - basegid + 1)){
48 std::cout << "Error: Map is invalid for test - fix test" << std::endl;
49 std::cerr << "Error: Map is invalid for test - fix test" << std::endl;
50 std::cout << "FAIL" << std::endl;
51 exit(1);
52 }
53 RCP<Array<zgno_t> > mygids = rcp(new Array<zgno_t>);
54 zgno_t firstzgno_t = proc;
55 if (firstzgno_t < basegid){
56 zgno_t n = basegid % proc;
57 if (n>0)
58 firstzgno_t = basegid - n + proc;
59 else
60 firstzgno_t = basegid;
61 }
62 for (zgno_t gid=firstzgno_t; gid <= maxgid; gid+=nprocs){
63 (*mygids).append(gid);
64 }
65
66 ArrayRCP<zgno_t> newIdArcp = Teuchos::arcp(mygids);
67
68 return newIdArcp;
69}
70
71ArrayRCP<zgno_t> roundRobinMap(const Tpetra::Map<zlno_t, zgno_t, znode_t> &tmap)
72{
73 const RCP<const Comm<int> > &comm = tmap.getComm();
74 int proc = comm->getRank();
75 int nprocs = comm->getSize();
76 zgno_t basegid = tmap.getMinAllGlobalIndex();
77 zgno_t maxgid = tmap.getMaxAllGlobalIndex();
78 size_t nglobalrows = tmap.getGlobalNumElements();
79
80 return roundRobinMapShared(proc, nprocs, basegid, maxgid, nglobalrows);
81}
82
83ArrayRCP<zgno_t> roundRobinMap(const Xpetra::Map<zlno_t, zgno_t, znode_t> &xmap)
84{
85 const RCP<const Comm<int> > &comm = xmap.getComm();
86 int proc = comm->getRank();
87 int nprocs = comm->getSize();
88 zgno_t basegid = xmap.getMinAllGlobalIndex();
89 zgno_t maxgid = xmap.getMaxAllGlobalIndex();
90 size_t nglobalrows = xmap.getGlobalNumElements();
91
92 return roundRobinMapShared(proc, nprocs, basegid, maxgid, nglobalrows);
93}
94
95int main(int narg, char *arg[])
96{
97 Tpetra::ScopeGuard tscope(&narg, &arg);
98 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
99
100 int rank = comm->getRank();
101 bool aok = true;
102
103 Teuchos::RCP<Teuchos::FancyOStream> outStream =
104 Teuchos::VerboseObjectBase::getDefaultOStream();
105 Teuchos::EVerbosityLevel v=Teuchos::VERB_EXTREME;
106
107 typedef Tpetra::CrsMatrix<zscalar_t,zlno_t,zgno_t,znode_t> tmatrix_t;
108 typedef Tpetra::CrsGraph<zlno_t,zgno_t,znode_t> tgraph_t;
109 typedef Tpetra::Vector<zscalar_t,zlno_t,zgno_t,znode_t> tvector_t;
110 typedef Tpetra::MultiVector<zscalar_t,zlno_t,zgno_t,znode_t> tmvector_t;
111 typedef Xpetra::CrsMatrix<zscalar_t,zlno_t,zgno_t,znode_t> xmatrix_t;
112 typedef Xpetra::CrsGraph<zlno_t,zgno_t,znode_t> xgraph_t;
113 typedef Xpetra::Vector<zscalar_t,zlno_t,zgno_t,znode_t> xvector_t;
114 typedef Xpetra::MultiVector<zscalar_t,zlno_t,zgno_t,znode_t> xmvector_t;
115
116 // Create object that can give us test Tpetra and Xpetra input.
117
118 RCP<UserInputForTests> uinput;
119 Teuchos::ParameterList params;
120 params.set("input file", "simple");
121 params.set("file type", "Chaco");
122
123 try{
124 uinput = rcp(new UserInputForTests(params, comm));
125 }
126 catch(std::exception &e){
127 aok = false;
128 std::cout << e.what() << std::endl;
129 }
130 TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
131
133 // Tpetra::CrsMatrix
134 // Tpetra::CrsGraph
135 // Tpetra::Vector
136 // Tpetra::MultiVector
138
139
140 // XpetraTraits<Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> >
141 {
142 RCP<tmatrix_t> M;
143
144 try{
145 M = uinput->getUITpetraCrsMatrix();
146 }
147 catch(std::exception &e){
148 aok = false;
149 std::cout << e.what() << std::endl;
150 }
151 TEST_FAIL_AND_EXIT(*comm, aok, "getTpetraCrsMatrix ", 1);
152
153 if (rank== 0)
154 std::cout << "Original Tpetra matrix " << M->getGlobalNumRows()
155 << " x " << M->getGlobalNumCols() << std::endl;
156
157 M->describe(*outStream,v);
158
159 ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(M->getRowMap()));
160
161 zgno_t localNumRows = newRowIds.size();
162
163 RCP<const tmatrix_t> newM;
164 try{
166 localNumRows, newRowIds.getRawPtr());
167 }
168 catch(std::exception &e){
169 aok = false;
170 std::cout << e.what() << std::endl;
171 }
172 TEST_FAIL_AND_EXIT(*comm, aok,
173 " Zoltan2::XpetraTraits<tmatrix_t>::doMigration ", 1);
174
175 if (rank== 0)
176 std::cout << "Migrated Tpetra matrix" << std::endl;
177
178 newM->describe(*outStream,v);
179 }
180
181 // XpetraTraits<Tpetra::CrsGraph<zscalar_t, zlno_t, zgno_t, znode_t> >
182 {
183 RCP<tgraph_t> G;
184
185 try{
186 G = uinput->getUITpetraCrsGraph();
187 }
188 catch(std::exception &e){
189 aok = false;
190 std::cout << e.what() << std::endl;
191 }
192 TEST_FAIL_AND_EXIT(*comm, aok, "getTpetraCrsGraph ", 1);
193
194 if (rank== 0)
195 std::cout << "Original Tpetra graph" << std::endl;
196
197 G->describe(*outStream,v);
198
199 ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(G->getRowMap()));
200
201 zgno_t localNumRows = newRowIds.size();
202
203 RCP<const tgraph_t> newG;
204 try{
206 localNumRows, newRowIds.getRawPtr());
207 }
208 catch(std::exception &e){
209 aok = false;
210 std::cout << e.what() << std::endl;
211 }
212 TEST_FAIL_AND_EXIT(*comm, aok,
213 " Zoltan2::XpetraTraits<tgraph_t>::doMigration ", 1);
214
215 if (rank== 0)
216 std::cout << "Migrated Tpetra graph" << std::endl;
217
218 newG->describe(*outStream,v);
219 }
220
221 // XpetraTraits<Tpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>>
222 {
223 RCP<tvector_t> V;
224
225 try{
226 V = rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 1));
227 V->randomize();
228 }
229 catch(std::exception &e){
230 aok = false;
231 std::cout << e.what() << std::endl;
232 }
233 TEST_FAIL_AND_EXIT(*comm, aok, "getTpetraVector", 1);
234
235 if (rank== 0)
236 std::cout << "Original Tpetra vector" << std::endl;
237
238 V->describe(*outStream,v);
239
240 ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(V->getMap()));
241
242 zgno_t localNumRows = newRowIds.size();
243
244 RCP<const tvector_t> newV;
245 try{
247 localNumRows, newRowIds.getRawPtr());
248 }
249 catch(std::exception &e){
250 aok = false;
251 std::cout << e.what() << std::endl;
252 }
253 TEST_FAIL_AND_EXIT(*comm, aok,
254 " Zoltan2::XpetraTraits<tvector_t>::doMigration ", 1);
255
256 if (rank== 0)
257 std::cout << "Migrated Tpetra vector" << std::endl;
258
259 newV->describe(*outStream,v);
260 }
261
262 // XpetraTraits<Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>>
263 {
264 RCP<tmvector_t> MV;
265
266 try{
267 MV = rcp(new tmvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 3));
268 MV->randomize();
269 }
270 catch(std::exception &e){
271 aok = false;
272 std::cout << e.what() << std::endl;
273 }
274 TEST_FAIL_AND_EXIT(*comm, aok, "getTpetraMultiVector", 1);
275
276 if (rank== 0)
277 std::cout << "Original Tpetra multivector" << std::endl;
278
279 MV->describe(*outStream,v);
280
281 ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(MV->getMap()));
282
283 zgno_t localNumRows = newRowIds.size();
284
285 RCP<const tmvector_t> newMV;
286 try{
288 localNumRows, newRowIds.getRawPtr());
289 }
290 catch(std::exception &e){
291 aok = false;
292 std::cout << e.what() << std::endl;
293 }
294 TEST_FAIL_AND_EXIT(*comm, aok,
295 " Zoltan2::XpetraTraits<tmvector_t>::doMigration ", 1);
296
297 if (rank== 0)
298 std::cout << "Migrated Tpetra multivector" << std::endl;
299
300 newMV->describe(*outStream,v);
301 }
302
304 // Xpetra::CrsMatrix
305 // Xpetra::CrsGraph
306 // Xpetra::Vector
307 // Xpetra::MultiVector
309
310 // XpetraTraits<Xpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> >
311 {
312 RCP<xmatrix_t> M;
313
314 try{
315 M = uinput->getUIXpetraCrsMatrix();
316 }
317 catch(std::exception &e){
318 aok = false;
319 std::cout << e.what() << std::endl;
320 }
321 TEST_FAIL_AND_EXIT(*comm, aok, "getXpetraCrsMatrix ", 1);
322
323 if (rank== 0)
324 std::cout << "Original Xpetra matrix" << std::endl;
325
326 M->describe(*outStream,v);
327
328 ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(M->getRowMap()));
329
330 zgno_t localNumRows = newRowIds.size();
331
332 RCP<const xmatrix_t> newM;
333 try{
335 localNumRows, newRowIds.getRawPtr());
336 }
337 catch(std::exception &e){
338 aok = false;
339 std::cout << e.what() << std::endl;
340 }
341 TEST_FAIL_AND_EXIT(*comm, aok,
342 " Zoltan2::XpetraTraits<xmatrix_t>::doMigration ", 1);
343
344 if (rank== 0)
345 std::cout << "Migrated Xpetra matrix" << std::endl;
346
347 newM->describe(*outStream,v);
348 }
349
350 // XpetraTraits<Xpetra::CrsGraph<zscalar_t, zlno_t, zgno_t, znode_t> >
351 {
352 RCP<xgraph_t> G;
353
354 try{
355 G = uinput->getUIXpetraCrsGraph();
356 }
357 catch(std::exception &e){
358 aok = false;
359 std::cout << e.what() << std::endl;
360 }
361 TEST_FAIL_AND_EXIT(*comm, aok, "getXpetraCrsGraph ", 1);
362
363 if (rank== 0)
364 std::cout << "Original Xpetra graph" << std::endl;
365
366 G->describe(*outStream,v);
367
368 ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(G->getRowMap()));
369
370 zgno_t localNumRows = newRowIds.size();
371
372 RCP<const xgraph_t> newG;
373 try{
375 localNumRows, newRowIds.getRawPtr());
376 }
377 catch(std::exception &e){
378 aok = false;
379 std::cout << e.what() << std::endl;
380 }
381 TEST_FAIL_AND_EXIT(*comm, aok,
382 " Zoltan2::XpetraTraits<xgraph_t>::doMigration ", 1);
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 aok = false;
402 std::cout << e.what() << std::endl;
403 }
404 TEST_FAIL_AND_EXIT(*comm, aok, "getXpetraVector", 1);
405
406 if (rank== 0)
407 std::cout << "Original Xpetra vector" << std::endl;
408
409 V->describe(*outStream,v);
410
411 ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(V->getMap()));
412
413 zgno_t localNumRows = newRowIds.size();
414
415 RCP<const xvector_t> newV;
416 try{
418 localNumRows, newRowIds.getRawPtr());
419 }
420 catch(std::exception &e){
421 aok = false;
422 std::cout << e.what() << std::endl;
423 }
424 TEST_FAIL_AND_EXIT(*comm, aok,
425 " Zoltan2::XpetraTraits<xvector_t>::doMigration ", 1);
426
427 if (rank== 0)
428 std::cout << "Migrated Xpetra vector" << std::endl;
429
430 newV->describe(*outStream,v);
431 }
432
433 // XpetraTraits<Xpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>>
434 {
435 RCP<xmvector_t> MV;
436
437 try{
438 RCP<tmvector_t> tMV =
439 rcp(new tmvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), 3));
440 tMV->randomize();
442 }
443 catch(std::exception &e){
444 aok = false;
445 std::cout << e.what() << std::endl;
446 }
447 TEST_FAIL_AND_EXIT(*comm, aok, "getXpetraMultiVector", 1);
448
449 if (rank== 0)
450 std::cout << "Original Xpetra multivector" << std::endl;
451
452 MV->describe(*outStream,v);
453
454 ArrayRCP<zgno_t> newRowIds = roundRobinMap(*(MV->getMap()));
455
456 zgno_t localNumRows = newRowIds.size();
457
458 RCP<const xmvector_t> newMV;
459 try{
461 localNumRows, newRowIds.getRawPtr());
462 }
463 catch(std::exception &e){
464 aok = false;
465 std::cout << e.what() << std::endl;
466 }
467 TEST_FAIL_AND_EXIT(*comm, aok,
468 " Zoltan2::XpetraTraits<xmvector_t>::doMigration ", 1);
469
470 if (rank== 0)
471 std::cout << "Migrated Xpetra multivector" << std::endl;
472
473 newMV->describe(*outStream,v);
474 }
475
477 // DONE
479
480 if (rank==0)
481 std::cout << "PASS" << std::endl;
482}
483
#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.