Zoltan2
Loading...
Searching...
No Matches
XpetraMultiVectorInput.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// @HEADER
11// ***********************************************************************
12// Zoltan2: Sandia Partitioning Ordering & Coloring Library
13// ***********************************************************************
14
20#include <string>
21
25
26#include <Teuchos_DefaultComm.hpp>
27#include <Teuchos_RCP.hpp>
28#include <Teuchos_Comm.hpp>
29#include <Teuchos_CommHelpers.hpp>
30
31using Teuchos::RCP;
32using Teuchos::rcp;
33using Teuchos::rcp_const_cast;
34using Teuchos::Comm;
35
36typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> tvector_t;
37typedef Xpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> xvector_t;
38
40template <typename User>
43{
44 RCP<const Comm<int> > comm = vector.getMap()->getComm();
45 int fail = 0, gfail=0;
46
47 if (!fail && ia.getNumEntriesPerID() !=nvec)
48 fail = 42;
49
50 size_t length = vector.getLocalLength();
51
52 if (!fail && ia.getLocalNumIDs() != length)
53 fail = 4;
54
55 gfail = globalFail(*comm, fail);
56
57 if (!gfail){
58 const zgno_t *vtxIds=NULL;
59 const zscalar_t *vals=NULL;
60 int stride;
61
62 size_t nvals = ia.getLocalNumIDs();
63 if (nvals != vector.getLocalLength())
64 fail = 8;
65
66 ia.getIDsView(vtxIds);
67
68 for (int v=0; v < nvec; v++){
69 auto vecdata = vector.getData(v);
70
71 ia.getEntriesView(vals, stride, v);
72
73 if (!fail && stride != 1)
74 fail = 10;
75
76 // check the values returned
77 for (size_t i = 0; i < length; i++)
78 if (vals[i*stride] != vecdata[i]) {
79 fail = 104;
80 }
81 }
82
83 gfail = globalFail(*comm, fail);
84 }
85
86
87 return gfail;
88}
89
90
91template <typename User>
94 std::vector<const zscalar_t *> &weights, std::vector<int> &strides)
95{
96 // Check the input adapter
97 int gfail = verifyInputAdapter(ia, vector, nvec);
98 int fail = gfail;
99
100 RCP<const Comm<int> > comm = vector.getMap()->getComm();
101 int wdim = weights.size();
102
103 // Check the input adapter weights
104 if (!gfail && (ia.getNumWeightsPerID() != wdim)) fail = 41;
105
106 gfail = globalFail(*comm, fail);
107
108 if (!gfail && wdim) {
109 const zscalar_t *wgt =NULL;
110 int stride;
111
112 for (int w=0; !fail && w < wdim; w++){
113 ia.getWeightsView(wgt, stride, w);
114
115 if (!fail && stride != strides[w])
116 fail = 101;
117
118 for (size_t v=0; !fail && v < vector.getLocalLength(); v++){
119 if (wgt[v*stride] != weights[w][v*stride])
120 fail=102;
121 }
122 }
123
124 gfail = globalFail(*comm, fail);
125 }
126 return gfail;
127}
128
129
131
133 std::ifstream &fp,
134 size_t &nIDs,
135 size_t &nEdges,
136 int &nWgts)
137{
138 // Read the header info from a Chaco .graph file
139 std::string line;
140 std::getline(fp, line);
141 while (line[0]=='#') std::getline(fp, line); // skip comments
142 std::istringstream issHeader(line);
143 std::string s_code;
144 issHeader >> nIDs >> nEdges >> s_code;
145 if (!strcmp(s_code.c_str(), "010") || !strcmp(s_code.c_str(), "011")) {
146 if (!(issHeader >> nWgts)) nWgts = 1;
147 }
148}
149
150template <typename User>
153 const char *fileprefixInp,
154 const Teuchos::Comm<int> &comm
155)
156{
157 // Compares files generated by generateFiles (in fileprefixGen.*)
158 // to input files (in fileprefixInp.*).
159
160 int fail = 0, gfail=0;
161
162 std::string tmp(fileprefixInp);
163 tmp = tmp + "_generated";
164 const char *fileprefixGen = tmp.c_str();
165
166 ia.generateFiles(fileprefixGen, comm);
167
168 // Only rank zero needs to check the resulting files
169 if (comm.getRank() == 0) {
170
171 size_t nIDsGen, nIDsInp;
172 size_t nEdgesGen, nEdgesInp;
173 int nWgtsGen = 0, nWgtsInp = 0;
174 std::string lineGen, lineInp;
175
176 std::ifstream fpGen, fpInp;
177 std::string graphFilenameGen = fileprefixGen;
178 graphFilenameGen = graphFilenameGen + ".graph";
179 std::string graphFilenameInp = fileprefixInp;
180 graphFilenameInp = graphFilenameInp + ".graph";
181
182 // Read header info from generated file
183 fpGen.open(graphFilenameGen.c_str(), std::ios::in);
184 readChacoGraphHeaderInfo(fpGen, nIDsGen, nEdgesGen, nWgtsGen);
185
186 // Read header info from input file
187 fpInp.open(graphFilenameInp.c_str(), std::ios::in);
188 readChacoGraphHeaderInfo(fpInp, nIDsInp, nEdgesInp, nWgtsInp);
189
190 // input file and generated file should have same number of IDs
191 if (nIDsGen != nIDsInp) {
192 std::cout << "GenerateFiles: nIDsGen " << nIDsGen
193 << " != nIDsInp " << nIDsInp << std::endl;
194 fail = 2222;
195 }
196
197 // Vector adapters don't have edges
198 if (!fail && nEdgesGen != 0) {
199 std::cout << "GenerateFiles: nEdgesGen " << nEdgesGen << " != 0"
200 << std::endl;
201 fail = 2223;
202 }
203
204 // Check the weights, if any
205 if (!fail && nWgtsGen) {
206 if (nWgtsInp) {
207 // input file has weights; compare weights
208 size_t cntWgtLines;
209 for (cntWgtLines = 0; cntWgtLines < nIDsGen &&
210 std::getline(fpGen, lineGen) &&
211 std::getline(fpInp, lineInp); cntWgtLines++) {
212
213 std::istringstream issGen(lineGen);
214 std::istringstream issInp(lineInp);
215
216 int nw = 0;
217 double wgtGen, wgtInp;
218
219 while (issGen >> wgtGen) {
220
221 if (nw < nWgtsInp) {
222 issInp >> wgtInp;
223 if (wgtGen != wgtInp) fail = 2224; // weights don't match
224 }
225 nw++;
226 }
227
228 if (nw != nWgtsGen) fail = 2225; // Too many or few weights on line
229 }
230 if (cntWgtLines != nIDsGen) fail = 2226; // not enough input lines
231 }
232 else {
233 // input file does not have weights;
234 // just make sure generated file has enough weights
235 size_t cntWgtLines = 0;
236 for (cntWgtLines = 0; cntWgtLines < nIDsGen &&
237 std::getline(fpGen, lineGen); cntWgtLines++) {
238 std::istringstream issGen(lineGen);
239 int nw = 0;
240 double wgtGen;
241 while (issGen) { issGen >> wgtGen; nw++; }
242 if (nw != nWgtsGen) fail = 2227; // Too many or few weights on line
243 }
244 if (cntWgtLines != nIDsGen) fail = 2228; // not enough input lines
245 }
246 }
247
248 fpGen.close();
249 fpInp.close();
250
251 // check coordinate files
252 if (!fail) {
253 std::string coordsFilenameGen = fileprefixGen;
254 coordsFilenameGen = coordsFilenameGen + ".coords";
255 std::string coordsFilenameInp = fileprefixInp;
256 coordsFilenameInp = coordsFilenameInp + ".coords";
257
258 fpGen.open(coordsFilenameGen.c_str(), std::ios::in);
259 fpInp.open(coordsFilenameInp.c_str(), std::ios::in);
260
261 size_t cnt;
262 for (cnt = 0; std::getline(fpGen,lineGen) &&
263 std::getline(fpInp,lineInp); cnt++) {
264
265 // Check each token
266 std::istringstream issGen(lineGen);
267 std::istringstream issInp(lineInp);
268
269 double xGen, xInp;
270 issGen >> xGen;
271 issInp >> xInp;
272 while (issGen && issInp) {
273
274 if (xGen != xInp) {
275 std::cout << "Coordinates " << xGen << " != " << xInp
276 << std::endl;
277 fail = 333;
278 }
279 issGen >> xGen;
280 issInp >> xInp;
281 }
282
283 // Check same number of tokens: are there any left in either line?
284 if (issGen || issInp) {
285 std::cout << "Dimension of generated file != dimension of input file"
286 << std::endl;
287 fail = 334;
288 }
289 }
290
291 // Did we have the correct number of coordinates?
292 if (!fail && cnt != nIDsGen) {
293 std::cout << "Number of coordinates read " << cnt
294 << " != number of IDs " << nIDsGen << std::endl;
295 fail = 444;
296 }
297
298 fpGen.close();
299 fpInp.close();
300 }
301 }
302
303 gfail = globalFail(comm, fail);
304 return gfail;
305}
306
308
309int main(int narg, char *arg[])
310{
311 Tpetra::ScopeGuard tscope(&narg, &arg);
312 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
313
314 int rank = comm->getRank();
315 int fail = 0, gfail=0;
316 bool aok = true;
317
318 // Read run-time options.
319 std::string inputFilePrefix("simple");
320 Teuchos::CommandLineProcessor cmdp (false, false);
321 cmdp.setOption("file", &inputFilePrefix,
322 "chaco file (prefix) to be used in test");
323 cmdp.parse(narg, arg);
324
325 // Create object that can give us test Tpetra, Xpetra
326 // and Epetra vectors for testing.
327
328 RCP<UserInputForTests> uinput;
329 Teuchos::ParameterList params;
330 params.set("input file", inputFilePrefix);
331 params.set("file type", "Chaco");
332
333 try{
334 uinput = rcp(new UserInputForTests(params, comm));
335 }
336 catch(std::exception &e){
337 aok = false;
338 std::cout << e.what() << std::endl;
339 }
340 TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
341
342 RCP<tvector_t> inputMVector; // original vector (for checking)
343 RCP<tvector_t> migratedMVector; // migrated vector
344
345 int nVec = 2;
346
347 inputMVector = uinput->getUICoordinates();
348 size_t vlen = inputMVector->getLocalLength();
349
350 // Vertex weights
351 std::vector<const zscalar_t *> IDWeights;
352 std::vector<int> IDWeightsStrides;
353
354 int nWeights = 0;
355 if (uinput->hasUIWeights()) {
356
357 auto uiweights = uinput->getUIWeights();
358
359 nWeights = uiweights->getNumVectors();
360 IDWeights.resize(nWeights);
361 IDWeightsStrides.resize(nWeights);
362
363 for (int w = 0; w < nWeights; w++) {
364 IDWeights[w] = uiweights->getData(w).getRawPtr();
365 IDWeightsStrides[w] = 1;
366 }
367 }
368
369 // To test migration in the input adapter we need a Solution object.
370
371 RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
372
375 typedef ia_t::part_t part_t;
376
377 part_t *p = new part_t [vlen];
378 memset(p, 0, sizeof(part_t) * vlen);
379 ArrayRCP<part_t> solnParts(p, 0, vlen, true);
380
381 soln_t solution(env, comm, nWeights);
382 solution.setParts(solnParts);
383
384
386 // User object is Tpetra::MultiVector, no weights
387 if (!gfail){
388 if (rank==0)
389 std::cout << "Constructed with Tpetra::MultiVector" << std::endl;
390
391 RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(inputMVector);
392 RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > tVInput;
393
394 try {
395 tVInput =
397 IDWeights, IDWeightsStrides));
398 }
399 catch (std::exception &e){
400 aok = false;
401 std::cout << e.what() << std::endl;
402 }
403 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter ", 1);
404
405 fail = verifyInputAdapter<tvector_t>(*tVInput, *inputMVector, nVec,
406 IDWeights, IDWeightsStrides);
407 fail = verifyGenerateFiles(*tVInput, inputFilePrefix.c_str(), *comm);
408
409 gfail = globalFail(*comm, fail);
410
411 if (!gfail){
412 tvector_t *vMigrate = NULL;
413 try{
414 tVInput->applyPartitioningSolution(*inputMVector, vMigrate, solution);
415 migratedMVector = rcp(vMigrate);
416 }
417 catch (std::exception &e){
418 fail = 11;
419 }
420
421 gfail = globalFail(*comm, fail);
422
423 if (!gfail){
424 RCP<const tvector_t> cnewV =
425 rcp_const_cast<const tvector_t>(migratedMVector);
426 RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
427 try{
428 newInput =
430 }
431 catch (std::exception &e){
432 aok = false;
433 std::cout << e.what() << std::endl;
434 }
435 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 2 ", 1);
436
437 if (rank==0){
438 std::cout << "Constructed with ";
439 std::cout << "Tpetra::MultiVector migrated to proc 0" << std::endl;
440 }
441 // We didn't apply partitioning soln to weights,
442 // so don't check them when verifying adapter
443 fail = verifyInputAdapter<tvector_t>(*newInput, *migratedMVector, nVec);
444 if (fail) fail += 100;
445 gfail = globalFail(*comm, fail);
446 }
447 }
448 if (gfail){
449 printFailureCode(*comm, fail);
450 }
451 }
452
454 // User object is Xpetra::MultiVector
455 if (!gfail){
456 if (rank==0)
457 std::cout << "Constructed with Xpetra::MultiVector" << std::endl;
458
459 RCP<xvector_t> xV =
461 RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
462 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
463
464 try {
465 xVInput =
467 IDWeights,
468 IDWeightsStrides));
469 }
470 catch (std::exception &e){
471 aok = false;
472 std::cout << e.what() << std::endl;
473 }
474 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 3 ", 1);
475
476 fail = verifyInputAdapter<xvector_t>(*xVInput, *inputMVector, nVec,
477 IDWeights, IDWeightsStrides);
478
479 gfail = globalFail(*comm, fail);
480
481 if (!gfail){
482 xvector_t *vMigrate =NULL;
483 try{
484 xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
485 }
486 catch (std::exception &e){
487 fail = 11;
488 }
489
490 gfail = globalFail(*comm, fail);
491
492 if (!gfail){
493 RCP<const xvector_t> cnewV(vMigrate);
494 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
495 try{
496 newInput =
498 }
499 catch (std::exception &e){
500 aok = false;
501 std::cout << e.what() << std::endl;
502 }
503 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 4 ", 1);
504
505 if (rank==0){
506 std::cout << "Constructed with ";
507 std::cout << "Xpetra::MultiVector migrated to proc 0" << std::endl;
508 }
509 fail = verifyInputAdapter<xvector_t>(*newInput, *migratedMVector, nVec);
510 if (fail) fail += 100;
511 gfail = globalFail(*comm, fail);
512 }
513 }
514 if (gfail){
515 printFailureCode(*comm, fail);
516 }
517 }
518
519#ifdef HAVE_EPETRA_DATA_TYPES
521 // User object is Epetra_MultiVector
522 typedef Epetra_MultiVector evector_t;
523 if (!gfail){
524 if (rank==0)
525 std::cout << "Constructed with Epetra_MultiVector" << std::endl;
526
527 RCP<evector_t> eV =
528 rcp(new Epetra_MultiVector(uinput->getUIEpetraCrsGraph()->RowMap(),
529 nVec));
530 for (int v = 0; v < nVec; v++) {
531 auto inV = inputMVector->getData(v);
532 for (int i = 0; i < eV->MyLength(); i++)
533 eV->ReplaceMyValue(i, v, inV[i]);
534 }
535
536 RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
537 RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
538
539 bool goodAdapter = true;
540 try {
541 eVInput =
543 IDWeights, IDWeightsStrides));
544 }
545 catch (std::exception &e){
546 if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
547 aok = false;
548 goodAdapter = false;
549 std::cout << e.what() << std::endl;
550 }
551 else {
552 // We expect an error from Xpetra when znode_t != Xpetra::EpetraNode
553 // Ignore it, but skip tests using this matrix adapter.
554 std::cout << "Node type is not supported by Xpetra's Epetra interface;"
555 << " Skipping this test." << std::endl;
556 std::cout << "FYI: Here's the exception message: " << std::endl
557 << e.what() << std::endl;
558 goodAdapter = false;
559 }
560 }
561 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 5 ", 1);
562
563 if (goodAdapter) {
564 fail = verifyInputAdapter<evector_t>(*eVInput, *inputMVector, nVec,
565 IDWeights, IDWeightsStrides);
566
567 gfail = globalFail(*comm, fail);
568
569 if (!gfail){
570 evector_t *vMigrate =NULL;
571 try{
572 eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
573 }
574 catch (std::exception &e){
575 fail = 11;
576 }
577
578 gfail = globalFail(*comm, fail);
579
580 if (!gfail){
581 RCP<const evector_t> cnewV(vMigrate, true);
582 RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
583 try{
584 newInput =
586 }
587 catch (std::exception &e){
588 aok = false;
589 std::cout << e.what() << std::endl;
590 }
591 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 6 ", 1);
592
593 if (rank==0){
594 std::cout << "Constructed with ";
595 std::cout << "Epetra_MultiVector migrated to proc 0" << std::endl;
596 }
597 fail = verifyInputAdapter<evector_t>(*newInput, *migratedMVector,
598 nVec);
599 if (fail) fail += 100;
600 gfail = globalFail(*comm, fail);
601 }
602 }
603 if (gfail){
604 printFailureCode(*comm, fail);
605 }
606 }
607 }
608#endif
609
611 // DONE
612
613 if (rank==0)
614 std::cout << "PASS" << std::endl;
615}
int globalFail(const Comm< int > &comm, int fail)
void printFailureCode(const Comm< int > &comm, int fail)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
int verifyInputAdapter(Zoltan2::XpetraMultiVectorAdapter< User > &ia, tvector_t &vector, int nvec)
Xpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
int verifyGenerateFiles(Zoltan2::VectorAdapter< User > &ia, const char *fileprefixInp, const Teuchos::Comm< int > &comm)
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
void readChacoGraphHeaderInfo(std::ifstream &fp, size_t &nIDs, size_t &nEdges, int &nWgts)
Traits for application input objects.
common code used by tests
float zscalar_t
Tpetra::Map ::global_ordinal_type zgno_t
Defines the XpetraMultiVectorAdapter.
int main()
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
A PartitioningSolution is a solution to a partitioning problem.
VectorAdapter defines the interface for vector input.
void generateFiles(const char *fileprefix, const Teuchos::Comm< int > &comm) const
Write files that can be used as input to Zoltan or Zoltan2 driver Creates chaco-formatted input files...
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
int getNumEntriesPerID() const
Return the number of vectors.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
size_t getLocalNumIDs() const
Returns the number of objects on this process.
static const std::string fail
static ArrayRCP< ArrayRCP< zscalar_t > > weights
SparseMatrixAdapter_t::part_t part_t
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.