Zoltan2
Loading...
Searching...
No Matches
UserInputForTests.hpp
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
14#ifndef USERINPUTFORTESTS
15#define USERINPUTFORTESTS
16
19#include <Zoltan2_Typedefs.hpp>
20
21#include <Tpetra_MultiVector.hpp>
22#include <Tpetra_CrsMatrix.hpp>
23#include <Tpetra_Map.hpp>
24#include <Xpetra_Vector.hpp>
25#include <Xpetra_CrsMatrix.hpp>
26#include <Xpetra_CrsGraph.hpp>
27
28#include <MatrixMarket_Tpetra.hpp>
29
30#ifdef HAVE_ZOLTAN2_GALERI
31#include <Galeri_XpetraProblemFactory.hpp>
32#include <Galeri_XpetraParameters.hpp>
33#endif
34
35#include <Tpetra_KokkosCompat_DefaultNode.hpp>
36
38#include <fstream>
39#include <string>
40
41#include <TpetraExt_MatrixMatrix_def.hpp>
42
43// pamgen required includes
45
46
47using Teuchos::RCP;
48using Teuchos::ArrayRCP;
49using Teuchos::ArrayView;
50using Teuchos::Array;
51using Teuchos::Comm;
52using Teuchos::rcp;
53using Teuchos::arcp;
54using Teuchos::rcp_const_cast;
55using Teuchos::ParameterList;
56using namespace Zoltan2_TestingFramework;
57
87
89{
90public:
91
92 typedef Tpetra::Map<zlno_t, zgno_t, znode_t> map_t;
93 typedef Tpetra::Export<zlno_t, zgno_t, znode_t> export_t;
94 typedef Tpetra::Import<zlno_t, zgno_t, znode_t> import_t;
95 typedef map_t::node_type default_znode_t;
96
97
117 UserInputForTests(string path, string testData,
118 const RCP<const Comm<int> > &c, bool debugInfo=false,
119 bool distributeInput=true);
120
137 UserInputForTests(int x, int y, int z, string matrixType,
138 const RCP<const Comm<int> > &c, bool debugInfo=false,
139 bool distributeInput=true);
140
156 UserInputForTests(const ParameterList &pList,
157 const RCP<const Comm<int> > &c);
158
161 static void getUIRandomData(unsigned int seed, zlno_t length,
162 zscalar_t min, zscalar_t max, ArrayView<ArrayRCP<zscalar_t > > data);
163
164 RCP<tMVector_t> getUICoordinates();
165
166 RCP<tMVector_t> getUIWeights();
167
168 RCP<tMVector_t> getUIEdgeWeights();
169
170 RCP<tcrsMatrix_t> getUITpetraCrsMatrix();
171
172 RCP<tcrsGraph_t> getUITpetraCrsGraph();
173
174 RCP<tVector_t> getUITpetraVector();
175
176 RCP<tMVector_t> getUITpetraMultiVector(int nvec);
177
178 RCP<xcrsMatrix_t> getUIXpetraCrsMatrix();
179
180 RCP<xcrsGraph_t> getUIXpetraCrsGraph();
181
182 RCP<xVector_t> getUIXpetraVector();
183
184 RCP<xMVector_t> getUIXpetraMultiVector(int nvec);
185
186#ifdef HAVE_ZOLTAN2_PAMGEN
187 PamgenMesh * getPamGenMesh(){return this->pamgen_mesh.operator->();}
188#endif
189
190 bool hasInput();
191
192 bool hasInputDataType(const string &input_type);
193
194 bool hasUICoordinates();
195
196 bool hasUIWeights();
197
198 bool hasUIEdgeWeights();
199
201
202 bool hasUITpetraCrsGraph();
203
204 bool hasUITpetraVector();
205
207
209
210 bool hasUIXpetraCrsGraph();
211
212 bool hasUIXpetraVector();
213
215
216 bool hasPamgenMesh();
217
218private:
219
220 bool verbose_;
221
222 const RCP<const Comm<int> > tcomm_;
223
224 bool havePamgenMesh;
225#ifdef HAVE_ZOLTAN2_PAMGEN
226 RCP<PamgenMesh> pamgen_mesh;
227#endif
228
229 RCP<tcrsMatrix_t> M_;
230 RCP<xcrsMatrix_t> xM_;
231
232 RCP<tMVector_t> xyz_;
233 RCP<tMVector_t> vtxWeights_;
234 RCP<tMVector_t> edgWeights_;
235
236 // Read a Matrix Market file into M_
237 // using Tpetra::MatrixMarket::Reader.
238 // If there are "Tim Davis" style coordinates
239 // that go with the file, read those into xyz_.
240
241 void readMatrixMarketFile(string path, string testData,bool distributeInput = true);
242
243 // Build matrix M_ from a mesh and a problem type
244 // with Galeri::Xpetra.
245
246 void buildCrsMatrix(int xdim, int ydim, int zdim, string type,
247 bool distributeInput);
248
249 // Read a Zoltan1 Chaco or Matrix Market file
250 // into M_. If it has geometric coordinates,
251 // read them into xyz_. If it has weights,
252 // read those into vtxWeights_ and edgWeights_.
253 void readZoltanTestData(string path, string testData,
254 bool distributeInput);
255
256 // Modify the Maps of an input matrix to make them non-contiguous
257 RCP<tcrsMatrix_t> modifyMatrixGIDs(RCP<tcrsMatrix_t> &in);
258 inline zgno_t newID(const zgno_t id) { return id * 2 + 10001; }
259
260 // Read Zoltan data that is in a .graph file.
261 void getUIChacoGraph(FILE *fptr, bool haveAssign, FILE *assignFile,
262 string name, bool distributeInput);
263
264 // Read Zoltan data that is in a .coords file.
265 void getUIChacoCoords(FILE *fptr, string name);
266
267 // Chaco reader code: This code is copied from zoltan/ch.
268 // It might benefit from a rewrite and simplification.
269
270 // Chaco reader helper functions: copied from zoltan/ch
271 static const int CHACO_LINE_LENGTH=200;
272 char chaco_line[CHACO_LINE_LENGTH]; /* space to hold values */
273 int chaco_offset; /* offset into line for next data */
274 int chaco_break_pnt; /* place in sequence to pause */
275 int chaco_save_pnt; /* place in sequence to save */
276
277 double chaco_read_val(FILE* infile, int *end_flag);
278 int chaco_read_int(FILE* infile, int *end_flag);
279 void chaco_flush_line(FILE*);
280
281 // Chaco graph reader: copied from zoltan/ch
282 int chaco_input_graph(FILE *fin, const char *inname, int **start,
283 int **adjacency, int *nvtxs, int *nVwgts,
284 float **vweights, int *nEwgts, float **eweights);
285
286 // Chaco coordinate reader: copied from zoltan/ch
287 int chaco_input_geom(FILE *fingeom, const char *geomname, int nvtxs,
288 int *igeom, double **x, double **y, double **z);
289
290 // Chaco coordinate reader: copied from zoltan/ch
291 int chaco_input_assign(FILE *finassign, const char *assignname, int nvtxs,
292 short *assignments);
293
294
295 // Read a GeomGen.txt file into M_
296 // Read coordinates into xyz_.
297 // If iti has weights read those to vtxWeights_
298 // and edgeWeights_
299 void readGeometricGenTestData(string path, string testData);
300
301 // Geometry Gnearatory helper function
302 void readGeoGenParams(string paramFileName,
303 ParameterList &geoparams);
304
305 // utility methods used when reading geom gen files
306
307 static string trim_right_copy(const string& s,
308 const string& delimiters = " \f\n\r\t\v" );
309
310 static string trim_left_copy(const string& s,
311 const string& delimiters = " \f\n\r\t\v" );
312
313 static string trim_copy(const string& s,
314 const string& delimiters = " \f\n\r\t\v" );
315
316
317 // Read a pamgen mesh
318 void readPamgenMeshFile(string path, string testData);
319#ifdef HAVE_ZOLTAN2_PAMGEN
320 void setPamgenAdjacencyGraph();
321 void setPamgenCoordinateMV();
322#endif
323};
324
325UserInputForTests::UserInputForTests(string path, string testData,
326 const RCP<const Comm<int> > &c,
327 bool debugInfo, bool distributeInput):
328verbose_(debugInfo), tcomm_(c), havePamgenMesh(false),
329M_(), xM_(), xyz_(), vtxWeights_(), edgWeights_(),
330chaco_offset(0), chaco_break_pnt(CHACO_LINE_LENGTH)
331{
332 bool zoltan1 = false;
333 string::size_type loc = path.find("/zoltan/test/"); // Zoltan1 data
334 if (loc != string::npos)
335 zoltan1 = true;
336
337 if (zoltan1)
338 readZoltanTestData(path, testData, distributeInput);
339 else
340 readMatrixMarketFile(path, testData);
341}
342
344 string matrixType,
345 const RCP<const Comm<int> > &c,
346 bool debugInfo,
347 bool distributeInput):
348verbose_(debugInfo), tcomm_(c), havePamgenMesh(false),
349M_(), xM_(), xyz_(), vtxWeights_(), edgWeights_(),
350chaco_offset(0), chaco_break_pnt(CHACO_LINE_LENGTH)
351{
352 if (matrixType.size() == 0){
353 int dim = 0;
354 if (x > 0) dim++;
355 if (y > 0) dim++;
356 if (z > 0) dim++;
357 if (dim == 1)
358 matrixType = string("Laplace1D");
359 else if (dim == 2)
360 matrixType = string("Laplace2D");
361 else if (dim == 3)
362 matrixType = string("Laplace3D");
363 else
364 throw std::runtime_error("input");
365
366 if (verbose_ && tcomm_->getRank() == 0)
367 std::cout << "UserInputForTests, Matrix type : " << matrixType << std::endl;
368 }
369
370 buildCrsMatrix(x, y, z, matrixType, distributeInput);
371}
372
373UserInputForTests::UserInputForTests(const ParameterList &pList,
374 const RCP<const Comm<int> > &c):
375tcomm_(c), havePamgenMesh(false),
376M_(), xM_(), xyz_(), vtxWeights_(), edgWeights_(),
377chaco_offset(0), chaco_break_pnt(CHACO_LINE_LENGTH)
378{
379
380 // get options
381 bool distributeInput = true, debugInfo = true;
382
383 if(pList.isParameter("distribute input"))
384 distributeInput = pList.get<bool>("distribute input");
385
386 if(pList.isParameter("debug"))
387 debugInfo = pList.get<bool>("debug");
388 this->verbose_ = debugInfo;
389
390 if(pList.isParameter("input file"))
391 {
392
393 // get input path
394 string path(".");
395 if(pList.isParameter("input path"))
396 path = pList.get<string>("input path");
397
398 string testData = pList.get<string>("input file");
399
400 // find out if we are working from the zoltan1 test diretory
402
403 // find out if we are using the geometric generator
404 if(pList.isParameter("file type") && pList.get<string>("file type") == "Geometric Generator")
405 file_format = GEOMGEN;
406 else if(pList.isParameter("file type") && pList.get<string>("file type") == "Pamgen")
407 {
408 file_format = PAMGEN;
409 }
410 else if(pList.isParameter("file type") && pList.get<string>("file type") == "Chaco")
411 file_format = CHACO; // this flag calls read ZoltanTestData, which calls the chaco readers...
412
413 // read the input file
414 switch (file_format) {
415 case GEOMGEN: readGeometricGenTestData(path,testData); break;
416 case PAMGEN: readPamgenMeshFile(path,testData); break;
417 case CHACO: readZoltanTestData(path, testData, distributeInput); break;
418 default: readMatrixMarketFile(path, testData, distributeInput); break;
419 }
420
421 }else if(pList.isParameter("x") || pList.isParameter("y") || pList.isParameter("z")){
422
423 int x,y,z;
424 x = y = z = 0;
425 if(pList.isParameter("x")) x = pList.get<int>("x");
426 if(pList.isParameter("y")) y = pList.get<int>("y");
427 if(pList.isParameter("z")) z = pList.get<int>("z");
428
429 string problemType = "";
430 if(pList.isParameter("equation type")) problemType = pList.get<string>("equation type");
431
432 if (problemType.size() == 0){
433 int dim = 0;
434 if (x > 0) dim++;
435 if (y > 0) dim++;
436 if (z > 0) dim++;
437 if (dim == 1)
438 problemType = string("Laplace1D");
439 else if (dim == 2)
440 problemType = string("Laplace2D");
441 else if (dim == 3)
442 problemType = string("Laplace3D");
443 else
444 throw std::runtime_error("input");
445
446 if (verbose_ && tcomm_->getRank() == 0)
447 std::cout << "UserInputForTests, Matrix type : " << problemType << std::endl;
448 }
449
450
451 buildCrsMatrix(x, y, z, problemType, distributeInput);
452
453 }else{
454 std::cerr << "Input file block undefined!" << std::endl;
455 }
456
457}
458
459
460RCP<Zoltan2_TestingFramework::tMVector_t> UserInputForTests::getUICoordinates()
461{
462 if (xyz_.is_null())
463 throw std::runtime_error("could not read coord file");
464 return xyz_;
465}
466
467RCP<Zoltan2_TestingFramework::tMVector_t> UserInputForTests::getUIWeights()
468{
469 return vtxWeights_;
470}
471
472RCP<Zoltan2_TestingFramework::tMVector_t> UserInputForTests::getUIEdgeWeights()
473{
474 return edgWeights_;
475}
476
477RCP<Zoltan2_TestingFramework::tcrsMatrix_t> UserInputForTests::getUITpetraCrsMatrix()
478{
479 if (M_.is_null())
480 throw std::runtime_error("could not read mtx file");
481 return M_;
482}
483
484RCP<Zoltan2_TestingFramework::tcrsGraph_t> UserInputForTests::getUITpetraCrsGraph()
485{
486 if (M_.is_null())
487 throw std::runtime_error("could not read mtx file");
488 return rcp_const_cast<tcrsGraph_t>(M_->getCrsGraph());
489}
490
491RCP<Zoltan2_TestingFramework::tVector_t> UserInputForTests::getUITpetraVector()
492{
493 RCP<tVector_t> V = rcp(new tVector_t(M_->getRowMap(), 1));
494 V->randomize();
495
496 return V;
497}
498
499RCP<Zoltan2_TestingFramework::tMVector_t> UserInputForTests::getUITpetraMultiVector(int nvec)
500{
501 RCP<tMVector_t> mV = rcp(new tMVector_t(M_->getRowMap(), nvec));
502 mV->randomize();
503
504 return mV;
505}
506
507RCP<Zoltan2_TestingFramework::xcrsMatrix_t> UserInputForTests::getUIXpetraCrsMatrix()
508{
509 if (M_.is_null())
510 throw std::runtime_error("could not read mtx file");
511 return xM_;
512}
513
514RCP<Zoltan2_TestingFramework::xcrsGraph_t> UserInputForTests::getUIXpetraCrsGraph()
515{
516 if (M_.is_null())
517 throw std::runtime_error("could not read mtx file");
518 return rcp_const_cast<xcrsGraph_t>(xM_->getCrsGraph());
519}
520
521RCP<Zoltan2_TestingFramework::xVector_t> UserInputForTests::getUIXpetraVector()
522{
524}
525
526RCP<Zoltan2_TestingFramework::xMVector_t> UserInputForTests::getUIXpetraMultiVector(int nvec)
527{
528 RCP<tMVector_t> tMV = getUITpetraMultiVector(nvec);
530}
531
533{
534 // find out if an input source has been loaded
535 return this->hasUICoordinates() || \
536 this->hasUITpetraCrsMatrix() || \
537 this->hasUITpetraCrsGraph() || \
538 this->hasPamgenMesh();
539}
540
541bool UserInputForTests::hasInputDataType(const string &input_type)
542{
543 if(input_type == "coordinates")
544 return this->hasUICoordinates();
545 else if(input_type == "tpetra_vector")
546 return this->hasUITpetraVector();
547 else if(input_type == "tpetra_multivector")
548 return this->hasUITpetraMultiVector();
549 else if(input_type == "tpetra_crs_graph")
550 return this->hasUITpetraCrsGraph();
551 else if(input_type == "tpetra_crs_matrix")
552 return this->hasUITpetraCrsMatrix();
553 else if(input_type == "xpetra_vector")
554 return this->hasUIXpetraVector();
555 else if(input_type == "xpetra_multivector")
556 return this->hasUIXpetraMultiVector();
557 else if(input_type == "xpetra_crs_graph")
558 return this->hasUIXpetraCrsGraph();
559 else if(input_type == "xpetra_crs_matrix")
560 return this->hasUIXpetraCrsMatrix();
561
562 return false;
563}
564
566{
567 return xyz_.is_null() ? false : true;
568}
569
571{
572 return vtxWeights_.is_null() ? false : true;
573}
574
576{
577 return edgWeights_.is_null() ? false : true;
578}
579
581{
582 return M_.is_null() ? false : true;
583}
584
586{
587 return M_.is_null() ? false : true;
588}
589
591{
592 return true;
593}
594
596{
597 return true;
598}
599
601{
602 return M_.is_null() ? false : true;
603}
604
606{
607 return M_.is_null() ? false : true;
608}
609
611{
612 return true;
613}
614
616{
617 return true;
618}
619
621{
622 return this->havePamgenMesh;
623}
624
625void UserInputForTests::getUIRandomData(unsigned int seed, zlno_t length,
626 zscalar_t min, zscalar_t max,
627 ArrayView<ArrayRCP<zscalar_t > > data)
628{
629 if (length < 1)
630 return;
631
632 size_t dim = data.size();
633 for (size_t i=0; i < dim; i++){
634 zscalar_t *tmp = new zscalar_t [length];
635 if (!tmp)
636 throw (std::bad_alloc());
637 data[i] = Teuchos::arcp(tmp, 0, length, true);
638 }
639
640 zscalar_t scalingFactor = (max-min) / RAND_MAX;
641 srand(seed);
642 for (size_t i=0; i < dim; i++){
643 zscalar_t *x = data[i].getRawPtr();
644 for (zlno_t j=0; j < length; j++)
645 *x++ = min + (zscalar_t(rand()) * scalingFactor);
646 }
647}
648
649// utility methods used when reading geom gen files
650
651string UserInputForTests::trim_right_copy(
652 const string& s,
653 const string& delimiters)
654{
655 return s.substr( 0, s.find_last_not_of( delimiters ) + 1 );
656}
657
658string UserInputForTests::trim_left_copy(
659 const string& s,
660 const string& delimiters)
661{
662 return s.substr( s.find_first_not_of( delimiters ) );
663}
664
665string UserInputForTests::trim_copy(
666 const string& s,
667 const string& delimiters)
668{
669 return trim_left_copy( trim_right_copy( s, delimiters ), delimiters );
670}
671
672void UserInputForTests::readGeometricGenTestData(string path,
673 string testData)
674{
675
676 std::ostringstream fname;
677 fname << path << "/" << testData << ".txt";
678
679 if (verbose_ && tcomm_->getRank() == 0)
680 std::cout << "UserInputForTests, Read: " << fname.str() << std::endl;
681
682 Teuchos::ParameterList geoparams("geo params");
683 readGeoGenParams(fname.str(),geoparams);
684
685 geometricgen_t * gg = new geometricgen_t(geoparams, this->tcomm_);
686
687 // get coordinate and point info
688 int coord_dim = gg->getCoordinateDimension();
689 int numWeightsPerCoord = gg->getNumWeights();
690 zlno_t numLocalPoints = gg->getNumLocalCoords();
691 zgno_t numGlobalPoints = gg->getNumGlobalCoords();
692
693 // allocate an array of coordinate arrays
694 zscalar_t **coords = new zscalar_t * [coord_dim];
695 for(int i = 0; i < coord_dim; ++i){
696 coords[i] = new zscalar_t[numLocalPoints];
697 }
698
699 // get a copy of the data
700 gg->getLocalCoordinatesCopy(coords);
701
702 // get an array of arrays of weight data (if any)
703 zscalar_t **weight = NULL;
704 if (numWeightsPerCoord) {
705 // memory allocation
706 weight = new zscalar_t * [numWeightsPerCoord];
707 for(int i = 0; i < numWeightsPerCoord; ++i){
708 weight[i] = new zscalar_t[numLocalPoints];
709 }
710
711 // get a copy of the weight data
712 gg->getLocalWeightsCopy(weight);
713 }
714
715 delete gg; // free up memory from geom gen
716
717
718 // make a Tpetra map
719 RCP<Tpetra::Map<zlno_t, zgno_t, znode_t> > mp =
720 rcp(new Tpetra::Map<zlno_t, zgno_t, znode_t>(numGlobalPoints, numLocalPoints, 0, this->tcomm_));
721
722 // make an array of array views containing the coordinate data
723 Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(coord_dim);
724 for (int i=0; i < coord_dim; i++){
725 if(numLocalPoints > 0){
726 Teuchos::ArrayView<const zscalar_t> a(coords[i], numLocalPoints);
727 coordView[i] = a;
728 }
729 else {
730 Teuchos::ArrayView<const zscalar_t> a;
731 coordView[i] = a;
732 }
733 }
734
735 // set the xyz_ multivector
736 xyz_ = RCP<tMVector_t>(new
737 tMVector_t(mp, coordView.view(0, coord_dim),
738 coord_dim));
739
740 // set the vtx weights
741 if (numWeightsPerCoord) {
742 // make an array of array views containing the weight data
743 Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > weightView(numWeightsPerCoord);
744 for (int i=0; i < numWeightsPerCoord; i++){
745 if(numLocalPoints > 0){
746 Teuchos::ArrayView<const zscalar_t> a(weight[i], numLocalPoints);
747 weightView[i] = a;
748 }
749 else {
750 Teuchos::ArrayView<const zscalar_t> a;
751 weightView[i] = a;
752 }
753 }
754
755 vtxWeights_ = RCP<tMVector_t>(new tMVector_t(mp, weightView.view(0, numWeightsPerCoord),
756 numWeightsPerCoord));
757 }
758}
759
760void UserInputForTests::readGeoGenParams(string paramFileName,
761 ParameterList &geoparams){
762
763 const char param_comment = '#';
764
765 std::string input = "";
766 char inp[25000];
767 for(int i = 0; i < 25000; ++i){
768 inp[i] = 0;
769 }
770
771 bool fail = false;
772 if(this->tcomm_->getRank() == 0){
773
774 std::fstream inParam(paramFileName.c_str());
775 if (inParam.fail())
776 {
777 fail = true;
778 }
779 if(!fail)
780 {
781 std::string tmp = "";
782 getline (inParam,tmp);
783 while (!inParam.eof()){
784 if(tmp != ""){
785 tmp = trim_copy(tmp);
786 if(tmp != ""){
787 input += tmp + "\n";
788 }
789 }
790 getline (inParam,tmp);
791 }
792 inParam.close();
793 for (size_t i = 0; i < input.size(); ++i){
794 inp[i] = input[i];
795 }
796 }
797 }
798
799
800
801 int size = (int)input.size();
802 if(fail){
803 size = -1;
804 }
805 this->tcomm_->broadcast(0, sizeof(int), (char*) &size);
806 if(size == -1){
807 throw "File " + paramFileName + " cannot be opened.";
808 }
809 this->tcomm_->broadcast(0, size, inp);
810 std::istringstream inParam(inp);
811 string str;
812 getline (inParam,str);
813 while (!inParam.eof()){
814 if(str[0] != param_comment){
815 size_t pos = str.find('=');
816 if(pos == string::npos){
817 throw "Invalid Line:" + str + " in parameter file";
818 }
819 string paramname = trim_copy(str.substr(0,pos));
820 string paramvalue = trim_copy(str.substr(pos + 1));
821 geoparams.set(paramname, paramvalue);
822 }
823 getline (inParam,str);
824 }
825}
826
828RCP<tcrsMatrix_t> UserInputForTests::modifyMatrixGIDs(
829 RCP<tcrsMatrix_t> &inMatrix
830)
831{
832 // Produce a new matrix with the same structure as inMatrix,
833 // but whose row/column GIDs are non-contiguous values.
834 // In this case GID g in inMatrix becomes g*2+1 in outMatrix.
835
836 // Create the map for the new matrix: same structure as inMap but with
837 // the GIDs modified.
838 RCP<const map_t> inMap = inMatrix->getRowMap();
839
840 size_t nRows = inMap->getLocalNumElements();
841 auto inRows = inMap->getMyGlobalIndices();
842 Teuchos::Array<zgno_t> outRows(nRows);
843 for (size_t i = 0; i < nRows; i++) {
844 outRows[i] = newID(inRows[i]);
845 }
846
847 Tpetra::global_size_t nGlobalRows = inMap->getGlobalNumElements();
848 RCP<map_t> outMap = rcp(new map_t(nGlobalRows, outRows(), 0,
849 inMap->getComm()));
850
851#ifdef INCLUDE_LENGTHY_OUTPUT
852 // Sanity check output
853 {
854 std::cout << inMap->getComm()->getRank() << " KDDKDD "
855 << "nGlobal " << inMap->getGlobalNumElements() << " "
856 << outMap->getGlobalNumElements() << "; "
857 << "nLocal " << inMap->getLocalNumElements() << " "
858 << outMap->getLocalNumElements() << "; "
859 << std::endl;
860 std::cout << inMap->getComm()->getRank() << " KDDKDD ";
861 for (size_t i = 0; i < nRows; i++)
862 std::cout << "(" << inMap->getMyGlobalIndices()[i] << ", "
863 << outMap->getMyGlobalIndices()[i] << ") ";
864 std::cout << std::endl;
865 }
866#endif // INCLUDE_LENGTHY_OUTPUT
867
868 // Create a new matrix using the new map
869 // Get the length of the longest row; allocate memory.
870 size_t rowLen = inMatrix->getLocalMaxNumRowEntries();
871 RCP<tcrsMatrix_t> outMatrix = rcp(new tcrsMatrix_t(outMap, rowLen));
872
873 typename tcrsMatrix_t::nonconst_global_inds_host_view_type indices("Indices", rowLen);
874 typename tcrsMatrix_t::nonconst_values_host_view_type values("Values", rowLen);
875
876 for (size_t i = 0; i < nRows; i++) {
877 size_t nEntries;
878 zgno_t inGid = inMap->getGlobalElement(i);
879 inMatrix->getGlobalRowCopy(inGid, indices, values, nEntries);
880 for (size_t j = 0; j < nEntries; j++)
881 indices[j] = newID(indices[j]);
882
883 zgno_t outGid = outMap->getGlobalElement(i);
884 ArrayView<zgno_t> Indices(indices.data(), nEntries);
885 ArrayView<zscalar_t> Values(values.data(), nEntries);
886 outMatrix->insertGlobalValues(outGid, Indices(0, nEntries),
887 Values(0, nEntries));
888 }
889 outMatrix->fillComplete();
890
891#ifdef INCLUDE_LENGTHY_OUTPUT
892 // Sanity check output
893 {
894 std::cout << inMap->getComm()->getRank() << " KDDKDD Rows "
895 << "nGlobal " << inMatrix->getGlobalNumRows() << " "
896 << outMatrix->getGlobalNumRows() << "; "
897 << "nLocal " << inMatrix->getLocalNumRows() << " "
898 << outMatrix->getLocalNumRows() << std::endl;
899 std::cout << inMap->getComm()->getRank() << " KDDKDD NNZS "
900 << "nGlobal " << inMatrix->getGlobalNumEntries() << " "
901 << outMatrix->getGlobalNumEntries() << "; "
902 << "nLocal " << inMatrix->getLocalNumEntries() << " "
903 << outMatrix->getLocalNumEntries() << std::endl;
904
905 size_t nIn, nOut;
906 Teuchos::Array<zgno_t> in(rowLen), out(rowLen);
907 Teuchos::Array<zscalar_t> inval(rowLen), outval(rowLen);
908
909 for (size_t i = 0; i < nRows; i++) {
910 std::cout << inMap->getComm()->getRank() << " KDDKDD " << i << " nnz(";
911 inMatrix->getGlobalRowCopy(inMap->getGlobalElement(i), in, inval, nIn);
912 outMatrix->getGlobalRowCopy(outMap->getGlobalElement(i), out, outval,
913 nOut);
914
915 std::cout << nIn << ", " << nOut << "): ";
916 for (size_t j = 0; j < nIn; j++) {
917 std::cout << "(" << in[j] << " " << inval[j] << ", "
918 << out[j] << " " << outval[j] << ") ";
919 }
920 std::cout << std::endl;
921 }
922 }
923#endif // INCLUDE_LENGTHY_OUTPUT
924
925 return outMatrix;
926}
927
929
930void UserInputForTests::readMatrixMarketFile(
931 string path,
932 string testData,
933 bool distributeInput
934)
935{
936 std::ostringstream fname;
937 fname << path << "/" << testData << ".mtx";
938
939 if (verbose_ && tcomm_->getRank() == 0)
940 std::cout << "UserInputForTests, Read: " << fname.str() << std::endl;
941
942 // FIXME (mfh 01 Aug 2016) Tpetra::MatrixMarket::Reader has a graph
943 // ("pattern" matrix) reader. Call its readSparseGraphFile method.
944
945 RCP<tcrsMatrix_t> toMatrix;
946 RCP<tcrsMatrix_t> fromMatrix;
947 bool aok = true;
948 try{
949 typedef Tpetra::MatrixMarket::Reader<tcrsMatrix_t> reader_type;
950 fromMatrix = reader_type::readSparseFile(fname.str(), tcomm_,
951 true, false, false);
952#ifdef KDD_NOT_READY_YET
953 // See note below about modifying coordinate IDs as well.
954 //if (makeNonContiguous)
955 fromMatrix = modifyMatrixGIDs(fromMatrix);
956#endif
957
958 if(!distributeInput)
959 {
960 if (verbose_ && tcomm_->getRank() == 0)
961 std::cout << "Constructing serial distribution of matrix" << std::endl;
962 // need to make a serial map and then import the data to redistribute it
963 RCP<const map_t> fromMap = fromMatrix->getRowMap();
964
965 size_t numGlobalCoords = fromMap->getGlobalNumElements();
966 size_t numLocalCoords = this->tcomm_->getRank() == 0 ? numGlobalCoords : 0;
967 RCP<const map_t> toMap = rcp(new map_t(numGlobalCoords,numLocalCoords, 0, tcomm_));
968
969 RCP<import_t> importer = rcp(new import_t(fromMap, toMap));
970 toMatrix = rcp(new tcrsMatrix_t(toMap,0));
971 toMatrix->doImport(*fromMatrix, *importer, Tpetra::INSERT);
972 toMatrix->fillComplete();
973
974 }else{
975 toMatrix = fromMatrix;
976 }
977 }catch (std::exception &e) {
978 if (tcomm_->getRank() == 0) {
979 std::cout << "UserInputForTests unable to read matrix market file:"
980 << fname.str() << std::endl;
981 std::cout << e.what() << std::endl;
982 }
983 aok = false;
984 }
985 TEST_FAIL_AND_THROW(*tcomm_, aok,
986 "UserInputForTests unable to read matrix market file");
987
988 M_ = toMatrix;
989#ifdef INCLUDE_LENGTHY_OUTPUT
990 std::cout << tcomm_->getRank() << " KDDKDD " << M_->getLocalNumRows()
991 << " " << M_->getGlobalNumRows()
992 << " " << M_->getLocalNumEntries()
993 << " " << M_->getGlobalNumEntries() << std::endl;
994#endif // INCLUDE_LENGTHY_OUTPUT
995
997
998 // Open the coordinate file.
999
1000 fname.str("");
1001 fname << path << "/" << testData << "_coord.mtx";
1002
1003 size_t coordDim = 0, numGlobalCoords = 0;
1004 size_t msg[2]={0,0};
1005 ArrayRCP<ArrayRCP<zscalar_t> > xyz;
1006 std::ifstream coordFile;
1007
1008 if (tcomm_->getRank() == 0){
1009
1010 if (verbose_)
1011 std::cout << "UserInputForTests, Read: " <<
1012 fname.str() << std::endl;
1013
1014 int fail = 0;
1015 try{
1016 coordFile.open(fname.str().c_str());
1017 }
1018 catch (std::exception &){ // there is no coordinate file
1019 fail = 1;
1020 }
1021
1022 if (!fail){
1023
1024 // Read past banner to number and dimension of coordinates.
1025
1026 char c[256];
1027 bool done=false;
1028
1029 while (!done && !fail && coordFile.good()){
1030 coordFile.getline(c, 256);
1031 if (!c[0])
1032 fail = 1;
1033 else if (c[0] == '%')
1034 continue;
1035 else {
1036 done=true;
1037 std::istringstream s(c);
1038 s >> numGlobalCoords >> coordDim;
1039 if (!s.eof() || numGlobalCoords < 1 || coordDim < 1)
1040 fail=1;
1041 }
1042 }
1043
1044 if (done){
1045
1046 // Read in the coordinates.
1047
1048 xyz = Teuchos::arcp(new ArrayRCP<zscalar_t> [coordDim], 0, coordDim);
1049
1050 for (size_t dim=0; !fail && dim < coordDim; dim++){
1051 size_t idx;
1052 zscalar_t *tmp = new zscalar_t [numGlobalCoords];
1053 if (!tmp)
1054 fail = 1;
1055 else{
1056 xyz[dim] = Teuchos::arcp(tmp, 0, numGlobalCoords);
1057
1058 for (idx=0; !coordFile.eof() && idx < numGlobalCoords; idx++){
1059 coordFile.getline(c, 256);
1060 std::istringstream s(c);
1061 s >> tmp[idx];
1062 }
1063
1064 if (idx < numGlobalCoords)
1065 fail = 1;
1066 }
1067 }
1068
1069 if (fail){
1070 ArrayRCP<zscalar_t> emptyArray;
1071 for (size_t dim=0; dim < coordDim; dim++)
1072 xyz[dim] = emptyArray; // free the memory
1073
1074 coordDim = 0;
1075 }
1076 }
1077 else{
1078 fail = 1;
1079 }
1080
1081 coordFile.close();
1082 }
1083
1084 msg[0] = coordDim;
1085 msg[1] = numGlobalCoords;
1086 }
1087
1088 // Broadcast coordinate dimension
1089 Teuchos::broadcast<int, size_t>(*tcomm_, 0, 2, msg);
1090
1091 coordDim = msg[0];
1092 numGlobalCoords = msg[1];
1093
1094 if (coordDim == 0)
1095 return;
1096
1097 zgno_t base = 0;
1098 RCP<const map_t> toMap;
1099
1100 if (!M_.is_null()){
1101 const RCP<const map_t> &mapM = M_->getRowMap();
1102 toMap = mapM;
1103 }
1104 else{
1105 if (verbose_ && tcomm_->getRank() == 0)
1106 {
1107 std::cout << "Matrix was null. ";
1108 std::cout << "Constructing distribution map for coordinate vector."
1109 << std::endl;
1110 }
1111
1112 if(!distributeInput)
1113 {
1114 if (verbose_ && tcomm_->getRank() == 0)
1115 std::cout << "Constructing serial distribution map for coordinates."
1116 << std::endl;
1117
1118 size_t numLocalCoords = this->tcomm_->getRank()==0 ? numGlobalCoords : 0;
1119 toMap = rcp(new map_t(numGlobalCoords,numLocalCoords, base, tcomm_));
1120 }else{
1121 toMap = rcp(new map_t(numGlobalCoords, base, tcomm_));
1122 }
1123 }
1124
1125 // Export coordinates to their owners
1126
1127 xyz_ = rcp(new tMVector_t(toMap, coordDim));
1128
1129 ArrayRCP<ArrayView<const zscalar_t> > coordLists(coordDim);
1130
1131 if (tcomm_->getRank() == 0){
1132
1133 for (size_t dim=0; dim < coordDim; dim++)
1134 coordLists[dim] = xyz[dim].view(0, numGlobalCoords);
1135
1136 zgno_t *tmp = new zgno_t [numGlobalCoords];
1137 if (!tmp)
1138 throw std::bad_alloc();
1139
1140 ArrayRCP<const zgno_t> rowIds = Teuchos::arcp(tmp, 0, numGlobalCoords);
1141
1142#ifdef KDD_NOT_READY_YET
1143 // TODO if modifyMatrixGIDs, we need to modify ids here as well
1144 for (zgno_t id=0; id < zgno_t(numGlobalCoords); id++)
1145 *tmp++ = newID(id);
1146#else
1147 for (zgno_t id=0; id < zgno_t(numGlobalCoords); id++)
1148 *tmp++ = id;
1149#endif
1150
1151 RCP<const map_t> fromMap = rcp(new map_t(numGlobalCoords,
1152 rowIds.view(0, numGlobalCoords),
1153 base, tcomm_));
1154
1155 tMVector_t allCoords(fromMap, coordLists.view(0, coordDim), coordDim);
1156
1157 export_t exporter(fromMap, toMap);
1158
1159 xyz_->doExport(allCoords, exporter, Tpetra::INSERT);
1160 }
1161 else{
1162
1163 RCP<const map_t> fromMap = rcp(new map_t(numGlobalCoords,
1164 ArrayView<zgno_t>(), base, tcomm_));
1165
1166 tMVector_t allCoords(fromMap, coordLists.view(0, coordDim), coordDim);
1167
1168 export_t exporter(fromMap, toMap);
1169
1170 xyz_->doExport(allCoords, exporter, Tpetra::INSERT);
1171 }
1172}
1173
1174void UserInputForTests::buildCrsMatrix(int xdim, int ydim, int zdim,
1175 string problemType, bool distributeInput)
1176{
1177#ifdef HAVE_ZOLTAN2_GALERI
1178 Teuchos::CommandLineProcessor tclp;
1179 Galeri::Xpetra::Parameters<zgno_t> params(tclp,
1180 xdim, ydim, zdim, problemType);
1181
1182 RCP<const Tpetra::Map<zlno_t, zgno_t> > map;
1183 if (distributeInput)
1184 map = rcp(new Tpetra::Map<zlno_t, zgno_t>(params.GetNumGlobalElements(),
1185 0, tcomm_));
1186 else {
1187 // All data initially on rank 0
1188 size_t nGlobalElements = params.GetNumGlobalElements();
1189 size_t nLocalElements = ((tcomm_->getRank() == 0) ? nGlobalElements : 0);
1190 map = rcp(new Tpetra::Map<zlno_t, zgno_t>(nGlobalElements, nLocalElements, 0,
1191 tcomm_));
1192 }
1193
1194 if (verbose_ && tcomm_->getRank() == 0){
1195
1196 std::cout << "Matrix is " << (distributeInput ? "" : "not");
1197 std::cout << "distributed." << std::endl;
1198
1199 std::cout << "UserInputForTests, Create matrix with " << problemType;
1200 std::cout << " (and " << xdim;
1201 if (zdim > 0)
1202 std::cout << " x " << ydim << " x " << zdim;
1203 else if (ydim > 0)
1204 std::cout << " x" << ydim << " x 1";
1205 else
1206 std::cout << "x 1 x 1";
1207
1208 std::cout << " mesh)" << std::endl;
1209
1210 }
1211
1212 bool aok = true;
1213 try{
1214 RCP<Galeri::Xpetra::Problem<Tpetra::Map<zlno_t, zgno_t>, Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t>, Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t> > > Pr =
1215 Galeri::Xpetra::BuildProblem<zscalar_t, zlno_t, zgno_t, Tpetra::Map<zlno_t, zgno_t>, Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t>, Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t> >
1216 (params.GetMatrixType(), map, params.GetParameterList());
1217 M_ = Pr->BuildMatrix();
1218 }
1219 catch (std::exception &) { // Probably not enough memory
1220 aok = false;
1221 }
1222 TEST_FAIL_AND_THROW(*tcomm_, aok,
1223 "UserInputForTests Galeri::Xpetra::BuildProblem failed");
1224
1226
1227 // Compute the coordinates for the matrix rows.
1228
1229 if (verbose_ && tcomm_->getRank() == 0)
1230 std::cout <<
1231 "UserInputForTests, Implied matrix row coordinates computed" <<
1232 std::endl;
1233
1234 ArrayView<const zgno_t> gids = map->getLocalElementList();
1235 zlno_t count = static_cast<zlno_t>(gids.size());
1236 int dim = 3;
1237 size_t pos = problemType.find("2D");
1238 if (pos != string::npos)
1239 dim = 2;
1240 else if (problemType == string("Laplace1D") ||
1241 problemType == string("Identity"))
1242 dim = 1;
1243
1244 Array<ArrayRCP<zscalar_t> > coordinates(dim);
1245
1246 if (count > 0){
1247 for (int i=0; i < dim; i++){
1248 zscalar_t *c = new zscalar_t [count];
1249 if (!c)
1250 throw(std::bad_alloc());
1251 coordinates[i] = Teuchos::arcp(c, 0, count, true);
1252 }
1253
1254 if (dim==3){
1255 zscalar_t *x = coordinates[0].getRawPtr();
1256 zscalar_t *y = coordinates[1].getRawPtr();
1257 zscalar_t *z = coordinates[2].getRawPtr();
1258 zgno_t xySize = xdim * ydim;
1259 for (zlno_t i=0; i < count; i++){
1260 zgno_t iz = gids[i] / xySize;
1261 zgno_t xy = gids[i] - iz*xySize;
1262 z[i] = zscalar_t(iz);
1263 y[i] = zscalar_t(xy / xdim);
1264 x[i] = zscalar_t(xy % xdim);
1265 }
1266 }
1267 else if (dim==2){
1268 zscalar_t *x = coordinates[0].getRawPtr();
1269 zscalar_t *y = coordinates[1].getRawPtr();
1270 for (zlno_t i=0; i < count; i++){
1271 y[i] = zscalar_t(gids[i] / xdim);
1272 x[i] = zscalar_t(gids[i] % xdim);
1273 }
1274 }
1275 else{
1276 zscalar_t *x = coordinates[0].getRawPtr();
1277 for (zlno_t i=0; i < count; i++)
1278 x[i] = zscalar_t(gids[i]);
1279 }
1280 }
1281
1282 Array<ArrayView<const zscalar_t> > coordView(dim);
1283 if (count > 0)
1284 for (int i=0; i < dim; i++)
1285 coordView[i] = coordinates[i].view(0,count);
1286
1287 xyz_ = rcp(new tMVector_t(map, coordView.view(0, dim), dim));
1288#else
1289 throw std::runtime_error("Galeri input requested but Trilinos is "
1290 "not built with Galeri.");
1291#endif
1292}
1293
1294void UserInputForTests::readZoltanTestData(string path, string testData,
1295 bool distributeInput)
1296{
1297 int rank = tcomm_->getRank();
1298 FILE *graphFile = NULL;
1299 FILE *coordFile = NULL;
1300 FILE *assignFile = NULL;
1301 int fileInfo[3];
1302
1303 for (int i = 0; i < CHACO_LINE_LENGTH; i++) chaco_line[i] = '\0';
1304
1305 if (rank == 0){
1306 // set chacho graph file name
1307 std::ostringstream chGraphFileName;
1308 chGraphFileName << path << "/" << testData << ".graph";
1309
1310 // set chaco graph
1311 std::ostringstream chCoordFileName;
1312 chCoordFileName << path << "/" << testData << ".coords";
1313
1314 // set chaco graph
1315 std::ostringstream chAssignFileName;
1316 chAssignFileName << path << "/" << testData << ".assign";
1317
1318 // open file
1319 graphFile = fopen(chGraphFileName.str().c_str(), "r");
1320
1321 if(!graphFile) // maybe the user is using the default zoltan1 path convention
1322 {
1323 chGraphFileName.str("");
1324 chCoordFileName.str("");
1325 // try constructing zoltan1 paths
1326 chGraphFileName << path << "/ch_" << testData << "/" << testData << ".graph";
1327 chCoordFileName << path << "/ch_" << testData << "/" << testData << ".coords";
1328 chAssignFileName << path << "/ch_" << testData << "/" << testData << ".assign";
1329 // try to open the graph file again, if this doesn't open
1330 // the user has not provided a valid path to the file
1331 graphFile = fopen(chGraphFileName.str().c_str(), "r");
1332 }
1333
1334 memset(fileInfo, 0, sizeof(int) * 3); // set fileinfo to 0's
1335 if (graphFile){
1336 fileInfo[0] = 1;
1337 if (verbose_ && tcomm_->getRank() == 0)
1338 std::cout << "UserInputForTests, open " <<
1339 chGraphFileName.str () << std::endl;
1340
1341 coordFile = fopen(chCoordFileName.str().c_str(), "r");
1342 if (coordFile){
1343 fileInfo[1] = 1;
1344 if (verbose_ && tcomm_->getRank() == 0)
1345 std::cout << "UserInputForTests, open " <<
1346 chCoordFileName.str () << std::endl;
1347 }
1348
1349 assignFile = fopen(chAssignFileName.str().c_str(), "r");
1350 if (assignFile){
1351 fileInfo[2] = 1;
1352 if (verbose_ && tcomm_->getRank() == 0)
1353 std::cout << "UserInputForTests, open " <<
1354 chAssignFileName.str () << std::endl;
1355 }
1356 }else{
1357 if (verbose_ && tcomm_->getRank() == 0){
1358 std::cout << "UserInputForTests, unable to open file: ";
1359 std::cout << chGraphFileName.str() << std::endl;
1360 }
1361 }
1362 }
1363
1364 // broadcast whether we have graphs and coords to all processes
1365 Teuchos::broadcast<int, int>(*tcomm_, 0, 3, fileInfo);
1366
1367 bool haveGraph = (fileInfo[0] == 1);
1368 bool haveCoords = (fileInfo[1] == 1);
1369 bool haveAssign = (fileInfo[2] == 1);
1370
1371 if (haveGraph){
1372 // builds M_, vtxWeights_, and edgWeights_ and closes file.
1373 try{
1374 getUIChacoGraph(graphFile, haveAssign, assignFile,
1375 testData, distributeInput);
1376 }
1378
1379 if (haveCoords){
1380 // builds xyz_ and closes the file.
1381 try{
1382 getUIChacoCoords(coordFile, testData);
1383 }
1385 }
1386 }
1387
1389}
1390
1391void UserInputForTests::getUIChacoGraph(FILE *fptr, bool haveAssign,
1392 FILE *assignFile, string fname,
1393 bool distributeInput)
1394{
1395 int rank = tcomm_->getRank();
1396 int graphCounts[5];
1397 int nvtxs=0, nedges=0;
1398 int nVwgts=0, nEwgts=0;
1399 int *start = NULL, *adj = NULL;
1400 float *ewgts = NULL, *vwgts = NULL;
1401 size_t *nzPerRow = NULL;
1402 size_t maxRowLen = 0;
1403 zgno_t base = 0;
1404 ArrayRCP<const size_t> rowSizes;
1405 int fail = 0;
1406 bool haveEdges = true;
1407
1408 if (rank == 0){
1409
1410 memset(graphCounts, 0, 5*sizeof(int));
1411
1412 // Reads in the file and closes it when done.
1413 fail = chaco_input_graph(fptr, fname.c_str(), &start, &adj,
1414 &nvtxs, &nVwgts, &vwgts, &nEwgts, &ewgts);
1415
1416 // There are Zoltan2 test graphs that have no edges.
1417
1418 // nEwgts must be 1 or 0 - add error
1419
1420 if (start == NULL)
1421 haveEdges = false;
1422
1423 if (verbose_)
1424 {
1425 std::cout << "UserInputForTests, " << nvtxs << " vertices,";
1426 if (haveEdges)
1427 std::cout << start[nvtxs] << " edges,";
1428 else
1429 std::cout << "no edges,";
1430 std::cout << nVwgts << " vertex weights, ";
1431 std::cout << nEwgts << " edge weights" << std::endl;
1432 }
1433
1434 if (nvtxs==0)
1435 fail = true;
1436
1437 if (fail){
1438 Teuchos::broadcast<int, int>(*tcomm_, 0, 5, graphCounts);
1439 throw std::runtime_error("Unable to read chaco file");
1440 }
1441
1442 if (haveEdges)
1443 nedges = start[nvtxs];
1444
1445 nzPerRow = new size_t [nvtxs];
1446 if (!nzPerRow)
1447 throw std::bad_alloc();
1448 rowSizes = arcp(nzPerRow, 0, nvtxs, true);
1449
1450 if (haveEdges){
1451 for (int i=0; i < nvtxs; i++){
1452 nzPerRow[i] = start[i+1] - start[i];
1453 if (nzPerRow[i] > maxRowLen)
1454 maxRowLen = nzPerRow[i];
1455 }
1456 }
1457 else{
1458 memset(nzPerRow, 0, sizeof(size_t) * nvtxs);
1459 }
1460
1461 // Make sure base gid is zero.
1462
1463 if (nedges){
1464 int chbase = 1; // chaco input files are one-based
1465 for (int i=0; i < nedges; i++)
1466 adj[i] -= chbase;
1467 }
1468
1469 graphCounts[0] = nvtxs;
1470 graphCounts[1] = nedges;
1471 graphCounts[2] = nVwgts;
1472 graphCounts[3] = nEwgts;
1473 graphCounts[4] = (int)maxRowLen; // size_t maxRowLen will fit; it is <= (int-int)
1474 }
1475
1476 Teuchos::broadcast<int, int>(*tcomm_, 0, 5, graphCounts);
1477
1478 if (graphCounts[0] == 0)
1479 throw std::runtime_error("Unable to read chaco file");
1480
1481 haveEdges = (graphCounts[1] > 0);
1482
1483 RCP<tcrsMatrix_t> fromMatrix;
1484 RCP<const map_t> fromMap;
1485
1486 // Create a Tpetra::CrsMatrix where rank 0 has entire matrix.
1487 if (rank == 0){
1488 fromMap = rcp(new map_t(nvtxs, nvtxs, base, tcomm_));
1489
1490 fromMatrix =
1491 rcp(new tcrsMatrix_t(fromMap, rowSizes()));
1492
1493 if (haveEdges){
1494
1495 zgno_t *edgeIds = new zgno_t [nedges];
1496 if (nedges && !edgeIds)
1497 throw std::bad_alloc();
1498 for (int i=0; i < nedges; i++)
1499 edgeIds[i] = adj[i];
1500
1501 free(adj);
1502 adj = NULL;
1503
1504 zgno_t *nextId = edgeIds;
1505 Array<zscalar_t> values(nedges, 1.0);
1506 if (nedges > 0 && nEwgts > 0) {
1507 for (int i=0; i < nedges; i++)
1508 values[i] = ewgts[i];
1509 free(ewgts);
1510 ewgts = NULL;
1511 }
1512
1513 for (int i=0; i < nvtxs; i++){
1514 if (nzPerRow[i] > 0){
1515 ArrayView<const zgno_t> rowNz(nextId, nzPerRow[i]);
1516 fromMatrix->insertGlobalValues(i, rowNz, values.view(start[i], start[i+1] - start[i]));
1517 nextId += nzPerRow[i];
1518 }
1519 }
1520
1521 delete [] edgeIds;
1522 edgeIds = NULL;
1523 }
1524
1525 fromMatrix->fillComplete();
1526 }
1527 else{
1528 nvtxs = graphCounts[0];
1529 nedges = graphCounts[1];
1530 nVwgts = graphCounts[2];
1531 nEwgts = graphCounts[3];
1532 maxRowLen = graphCounts[4];
1533
1534 // Create a Tpetra::CrsMatrix where rank 0 has entire matrix.
1535
1536 fromMap = rcp(new map_t(nvtxs, 0, base, tcomm_));
1537
1538 fromMatrix =
1539 rcp(new tcrsMatrix_t(fromMap, rowSizes()));
1540
1541 fromMatrix->fillComplete();
1542 }
1543
1544#ifdef KDDKDDPRINT
1545 if (rank == 0) {
1546 size_t sz = fromMatrix->getLocalMaxNumRowEntries();
1547 Teuchos::Array<zgno_t> indices(sz);
1548 Teuchos::Array<zscalar_t> values(sz);
1549 for (size_t i = 0; i < fromMatrix->getLocalNumRows(); i++) {
1550 zgno_t gid = fromMatrix->getRowMap()->getGlobalElement(i);
1551 size_t num;
1552 fromMatrix->getGlobalRowCopy(gid, indices(), values(), num);
1553 std::cout << "ROW " << gid << ": ";
1554 for (size_t j = 0; j < num; j++)
1555 std::cout << indices[j] << " ";
1556 std::cout << std::endl;
1557 }
1558 }
1559#endif
1560
1561 RCP<const map_t> toMap;
1562 RCP<tcrsMatrix_t> toMatrix;
1563 RCP<import_t> importer;
1564
1565 if (distributeInput) {
1566 if (haveAssign) {
1567 // Read assignments from Chaco assignment file
1568 short *assignments = new short[nvtxs];
1569 if (rank == 0) {
1570 fail = chaco_input_assign(assignFile, fname.c_str(), nvtxs, assignments);
1571 }
1572 // Broadcast coordinate dimension
1573 Teuchos::broadcast<int, short>(*tcomm_, 0, nvtxs, assignments);
1574
1575 // Build map with my vertices
1576 Teuchos::Array<zgno_t> mine;
1577 for (int i = 0; i < nvtxs; i++) {
1578 if (assignments[i] == rank)
1579 mine.push_back(i);
1580 }
1581
1582 Tpetra::global_size_t dummy =
1583 Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid();
1584 toMap = rcp(new map_t(dummy, mine(), base, tcomm_));
1585 delete [] assignments;
1586 }
1587 else {
1588 // Create a Tpetra::Map with default row distribution.
1589 toMap = rcp(new map_t(nvtxs, base, tcomm_));
1590 }
1591 toMatrix = rcp(new tcrsMatrix_t(toMap, maxRowLen));
1592
1593 // Import the data.
1594 importer = rcp(new import_t(fromMap, toMap));
1595 toMatrix->doImport(*fromMatrix, *importer, Tpetra::INSERT);
1596 toMatrix->fillComplete();
1597 }
1598 else {
1599 toMap = fromMap;
1600 toMatrix = fromMatrix;
1601 }
1602
1603 M_ = toMatrix;
1604
1605 // Vertex weights, if any
1606
1607 typedef ArrayRCP<const ArrayView<const zscalar_t> > arrayArray_t;
1608
1609 if (nVwgts > 0){
1610
1611 ArrayRCP<zscalar_t> weightBuf;
1612 ArrayView<const zscalar_t> *wgts = new ArrayView<const zscalar_t> [nVwgts];
1613
1614 if (rank == 0){
1615 size_t len = nVwgts * nvtxs;
1616 zscalar_t *buf = new zscalar_t [len];
1617 if (!buf) throw std::bad_alloc();
1618 weightBuf = arcp(buf, 0, len, true);
1619
1620 for (int widx=0; widx < nVwgts; widx++){
1621 wgts[widx] = ArrayView<const zscalar_t>(buf, nvtxs);
1622 float *vw = vwgts + widx;
1623 for (int i=0; i < nvtxs; i++, vw += nVwgts)
1624 buf[i] = *vw;
1625 buf += nvtxs;
1626 }
1627
1628 free(vwgts);
1629 vwgts = NULL;
1630 }
1631
1632 arrayArray_t vweights = arcp(wgts, 0, nVwgts, true);
1633
1634 RCP<tMVector_t> fromVertexWeights =
1635 rcp(new tMVector_t(fromMap, vweights.view(0, nVwgts), nVwgts));
1636
1637 RCP<tMVector_t> toVertexWeights;
1638 if (distributeInput) {
1639 toVertexWeights = rcp(new tMVector_t(toMap, nVwgts));
1640 toVertexWeights->doImport(*fromVertexWeights, *importer, Tpetra::INSERT);
1641 }
1642 else
1643 toVertexWeights = fromVertexWeights;
1644
1645 vtxWeights_ = toVertexWeights;
1646 }
1647
1648 // Edge weights, if any
1649
1650 if (haveEdges && nEwgts > 0){
1651
1652 // No longer distributing edge weights; they are the matrix values
1653 /*
1654 ArrayRCP<zscalar_t> weightBuf;
1655 ArrayView<const zscalar_t> *wgts = new ArrayView<const zscalar_t> [nEwgts];
1656
1657 toMap = rcp(new map_t(nedges, M_->getLocalNumEntries(), base, tcomm_));
1658
1659 if (rank == 0){
1660 size_t len = nEwgts * nedges;
1661 zscalar_t *buf = new zscalar_t [len];
1662 if (!buf) throw std::bad_alloc();
1663 weightBuf = arcp(buf, 0, len, true);
1664
1665 for (int widx=0; widx < nEwgts; widx++){
1666 wgts[widx] = ArrayView<const zscalar_t>(buf, nedges);
1667 float *ew = ewgts + widx;
1668 for (int i=0; i < nedges; i++, ew += nEwgts)
1669 buf[i] = *ew;
1670 buf += nedges;
1671 }
1672
1673 free(ewgts);
1674 ewgts = NULL;
1675 fromMap = rcp(new map_t(nedges, nedges, base, tcomm_));
1676 }
1677 else{
1678 fromMap = rcp(new map_t(nedges, 0, base, tcomm_));
1679 }
1680
1681 arrayArray_t eweights = arcp(wgts, 0, nEwgts, true);
1682
1683 RCP<tMVector_t> fromEdgeWeights;
1684 RCP<tMVector_t> toEdgeWeights;
1685 RCP<import_t> edgeImporter;
1686
1687 if (distributeInput) {
1688 fromEdgeWeights =
1689 rcp(new tMVector_t(fromMap, eweights.view(0, nEwgts), nEwgts));
1690 toEdgeWeights = rcp(new tMVector_t(toMap, nEwgts));
1691 edgeImporter = rcp(new import_t(fromMap, toMap));
1692 toEdgeWeights->doImport(*fromEdgeWeights, *edgeImporter, Tpetra::INSERT);
1693 }
1694 else
1695 toEdgeWeights = fromEdgeWeights;
1696
1697 edgWeights_ = toEdgeWeights;
1698 */
1699
1700 toMap = rcp(new map_t(nedges, M_->getLocalNumEntries(), base, tcomm_));
1701 edgWeights_ = rcp(new tMVector_t(toMap, nEwgts));
1702
1703 size_t maxSize = M_->getLocalMaxNumRowEntries();
1704 tcrsMatrix_t::nonconst_local_inds_host_view_type colind("colind", maxSize);
1705 tcrsMatrix_t::nonconst_values_host_view_type vals("values", maxSize);
1706 size_t nEntries;
1707
1708 for (size_t i = 0, idx = 0; i < M_->getLocalNumRows(); i++) {
1709 M_->getLocalRowCopy(i, colind, vals, nEntries);
1710 for (size_t j = 0; j < nEntries; j++) {
1711 edgWeights_->replaceLocalValue(idx, 0, vals[j]); // Assuming nEwgts==1
1712 idx++;
1713 }
1714 }
1715 }
1716
1717 if (start) {
1718 free(start);
1719 start = NULL;
1720 }
1721}
1722
1723void UserInputForTests::getUIChacoCoords(FILE *fptr, string fname)
1724{
1725 int rank = tcomm_->getRank();
1726 int ndim=0;
1727 double *x=NULL, *y=NULL, *z=NULL;
1728 int fail = 0;
1729
1730 size_t globalNumVtx = M_->getGlobalNumRows();
1731
1732 if (rank == 0){
1733
1734 // This function is from the Zoltan C-library.
1735
1736 // Reads in the file and closes it when done.
1737 fail = chaco_input_geom(fptr, fname.c_str(), (int)globalNumVtx,
1738 &ndim, &x, &y, &z);
1739
1740 if (fail)
1741 ndim = 0;
1742
1743 if (verbose_){
1744 std::cout << "UserInputForTests, read " << globalNumVtx;
1745 std::cout << " " << ndim << "-dimensional coordinates." << std::endl;
1746 }
1747 }
1748
1749 Teuchos::broadcast<int, int>(*tcomm_, 0, 1, &ndim);
1750
1751 if (ndim == 0)
1752 throw std::runtime_error("Can't read coordinate file");
1753
1754 ArrayRCP<ArrayRCP<const zscalar_t> > coords(ndim);
1755 zlno_t len = 0;
1756
1757 if (rank == 0){
1758
1759 for (int dim=0; dim < ndim; dim++){
1760 zscalar_t *v = new zscalar_t [globalNumVtx];
1761 if (!v)
1762 throw std::bad_alloc();
1763 coords[dim] = arcp<const zscalar_t>(v, 0, globalNumVtx, true);
1764 double *val = (dim==0 ? x : (dim==1 ? y : z));
1765 for (size_t i=0; i < globalNumVtx; i++)
1766 v[i] = zscalar_t(val[i]);
1767
1768 free(val);
1769 }
1770
1771 len = static_cast<zlno_t>(globalNumVtx);
1772 }
1773
1774 RCP<const map_t> fromMap = rcp(new map_t(globalNumVtx, len, 0, tcomm_));
1775 RCP<const map_t> toMap = M_->getRowMap();
1776 RCP<import_t> importer = rcp(new import_t(fromMap, toMap));
1777
1778 Array<ArrayView<const zscalar_t> > coordData;
1779 for (int dim=0; dim < ndim; dim++)
1780 coordData.push_back(coords[dim].view(0, len));
1781
1782 RCP<tMVector_t> fromCoords =
1783 rcp(new tMVector_t(fromMap, coordData.view(0, ndim), ndim));
1784
1785 RCP<tMVector_t> toCoords = rcp(new tMVector_t(toMap, ndim));
1786
1787 toCoords->doImport(*fromCoords, *importer, Tpetra::INSERT);
1788
1789 xyz_ = toCoords;
1790
1791}
1792
1795
1796double UserInputForTests::chaco_read_val(
1797 FILE* infile, /* file to read value from */
1798 int *end_flag /* 0 => OK, 1 => EOL, -1 => EOF */
1799)
1800{
1801 double val; /* return value */
1802 char *ptr; /* ptr to next string to read */
1803 char *ptr2; /* ptr to next string to read */
1804 int length; /* length of line to read */
1805 int length_left;/* length of line still around */
1806 int white_seen; /* have I detected white space yet? */
1807 int done; /* checking for end of scan */
1808 int i; /* loop counter */
1809
1810 *end_flag = 0;
1811
1812 if (chaco_offset == 0 || chaco_offset >= chaco_break_pnt) {
1813 if (chaco_offset >= chaco_break_pnt) { /* Copy rest of line back to beginning. */
1814 length_left = CHACO_LINE_LENGTH - chaco_save_pnt - 1;
1815 ptr2 = chaco_line;
1816 ptr = &chaco_line[chaco_save_pnt];
1817 for (i=length_left; i; i--) *ptr2++ = *ptr++;
1818 length = chaco_save_pnt + 1;
1819 }
1820 else {
1821 length = CHACO_LINE_LENGTH;
1822 length_left = 0;
1823 }
1824
1825 /* Now read next line, or next segment of current one. */
1826 ptr2 = fgets(&chaco_line[length_left], length, infile);
1827
1828 if (ptr2 == (char *) NULL) { /* We've hit end of file. */
1829 *end_flag = -1;
1830 return((double) 0.0);
1831 }
1832
1833 if ((chaco_line[CHACO_LINE_LENGTH - 2] != '\n') && (chaco_line[CHACO_LINE_LENGTH - 2] != '\f')
1834 && (strlen(chaco_line) == CHACO_LINE_LENGTH - 1)){
1835 /* Line too long. Find last safe place in chaco_line. */
1836 chaco_break_pnt = CHACO_LINE_LENGTH - 1;
1837 chaco_save_pnt = chaco_break_pnt;
1838 white_seen = 0;
1839 done = 0;
1840 while (!done) {
1841 --chaco_break_pnt;
1842 if (chaco_line[chaco_break_pnt] != '\0') {
1843 if (isspace((int)(chaco_line[chaco_break_pnt]))) {
1844 if (!white_seen) {
1845 chaco_save_pnt = chaco_break_pnt + 1;
1846 white_seen = 1;
1847 }
1848 }
1849 else if (white_seen) {
1850 done= 1;
1851 }
1852 }
1853 }
1854 }
1855 else {
1856 chaco_break_pnt = CHACO_LINE_LENGTH;
1857 }
1858
1859 chaco_offset = 0;
1860 }
1861
1862 while (isspace((int)(chaco_line[chaco_offset])) && chaco_offset < CHACO_LINE_LENGTH) chaco_offset++;
1863 if (chaco_line[chaco_offset] == '%' || chaco_line[chaco_offset] == '#') {
1864 *end_flag = 1;
1865 if (chaco_break_pnt < CHACO_LINE_LENGTH) {
1866 chaco_flush_line(infile);
1867 }
1868 return((double) 0.0);
1869 }
1870
1871 ptr = &(chaco_line[chaco_offset]);
1872 val = strtod(ptr, &ptr2);
1873
1874 if (ptr2 == ptr) { /* End of input line. */
1875 chaco_offset = 0;
1876 *end_flag = 1;
1877 return((double) 0.0);
1878 }
1879 else {
1880 chaco_offset = (int) (ptr2 - chaco_line) / sizeof(char);
1881 }
1882
1883 return(val);
1884}
1885
1886
1887int UserInputForTests::chaco_read_int(
1888 FILE *infile, /* file to read value from */
1889 int *end_flag /* 0 => OK, 1 => EOL, -1 => EOF */
1890)
1891{
1892 int val; /* return value */
1893 char *ptr; /* ptr to next string to read */
1894 char *ptr2; /* ptr to next string to read */
1895 int length; /* length of line to read */
1896 int length_left; /* length of line still around */
1897 int white_seen; /* have I detected white space yet? */
1898 int done; /* checking for end of scan */
1899 int i; /* loop counter */
1900
1901 *end_flag = 0;
1902
1903 if (chaco_offset == 0 || chaco_offset >= chaco_break_pnt) {
1904 if (chaco_offset >= chaco_break_pnt) { /* Copy rest of line back to beginning. */
1905 length_left = CHACO_LINE_LENGTH - chaco_save_pnt - 1;
1906 ptr2 = chaco_line;
1907 ptr = &chaco_line[chaco_save_pnt];
1908 for (i=length_left; i; i--) *ptr2++ = *ptr++;
1909 length = chaco_save_pnt + 1;
1910 }
1911 else {
1912 length = CHACO_LINE_LENGTH;
1913 length_left = 0;
1914 }
1915
1916 /* Now read next line, or next segment of current one. */
1917 ptr2 = fgets(&chaco_line[length_left], length, infile);
1918
1919 if (ptr2 == (char *) NULL) { /* We've hit end of file. */
1920 *end_flag = -1;
1921 return(0);
1922 }
1923
1924 if ((chaco_line[CHACO_LINE_LENGTH - 2] != '\n') && (chaco_line[CHACO_LINE_LENGTH - 2] != '\f')
1925 && (strlen(chaco_line) == CHACO_LINE_LENGTH - 1)){
1926 /* Line too long. Find last safe place in line. */
1927 chaco_break_pnt = CHACO_LINE_LENGTH - 1;
1928 chaco_save_pnt = chaco_break_pnt;
1929 white_seen = 0;
1930 done = 0;
1931 while (!done) {
1932 --chaco_break_pnt;
1933 if (chaco_line[chaco_break_pnt] != '\0') {
1934 if (isspace((int)(chaco_line[chaco_break_pnt]))) {
1935 if (!white_seen) {
1936 chaco_save_pnt = chaco_break_pnt + 1;
1937 white_seen = 1;
1938 }
1939 }
1940 else if (white_seen) {
1941 done= 1;
1942 }
1943 }
1944 }
1945 }
1946 else {
1947 chaco_break_pnt = CHACO_LINE_LENGTH;
1948 }
1949
1950 chaco_offset = 0;
1951 }
1952
1953 while (isspace((int)(chaco_line[chaco_offset])) && chaco_offset < CHACO_LINE_LENGTH) chaco_offset++;
1954 if (chaco_line[chaco_offset] == '%' || chaco_line[chaco_offset] == '#') {
1955 *end_flag = 1;
1956 if (chaco_break_pnt < CHACO_LINE_LENGTH) {
1957 chaco_flush_line(infile);
1958 }
1959 return(0);
1960 }
1961
1962 ptr = &(chaco_line[chaco_offset]);
1963 val = (int) strtol(ptr, &ptr2, 10);
1964
1965 if (ptr2 == ptr) { /* End of input chaco_line. */
1966 chaco_offset = 0;
1967 *end_flag = 1;
1968 return(0);
1969 }
1970 else {
1971 chaco_offset = (int) (ptr2 - chaco_line) / sizeof(char);
1972 }
1973
1974 return(val);
1975}
1976
1977void UserInputForTests::chaco_flush_line(
1978 FILE *infile /* file to read value from */
1979)
1980{
1981 char c; /* character being read */
1982
1983 c = fgetc(infile);
1984 while (c != '\n' && c != '\f')
1985 c = fgetc(infile);
1986}
1987
1988int UserInputForTests::chaco_input_graph(
1989 FILE *fin, /* input file */
1990 const char *inname, /* name of input file */
1991 int **start, /* start of edge list for each vertex */
1992 int **adjacency, /* edge list data */
1993 int *nvtxs, /* number of vertices in graph */
1994 int *nVwgts, /* # of vertex weights per node */
1995 float **vweights, /* vertex weight list data */
1996 int *nEwgts, /* # of edge weights per edge */
1997 float **eweights /* edge weight list data */
1998)
1999{
2000 int *adjptr; /* loops through adjacency data */
2001 float *ewptr; /* loops through edge weight data */
2002 int narcs; /* number of edges expected in graph */
2003 int nedges; /* twice number of edges really in graph */
2004 int nedge; /* loops through edges for each vertex */
2005 int flag; /* condition indicator */
2006 int skip_flag; /* should this edge be ignored? */
2007 int end_flag; /* indicates end of line or file */
2008 int vtx; /* vertex in graph */
2009 int line_num; /* line number in input file */
2010 int sum_edges; /* total number of edges read so far */
2011 int option = 0; /* input option */
2012 int using_ewgts; /* are edge weights in input file? */
2013 int using_vwgts; /* are vertex weights in input file? */
2014 int vtxnums; /* are vertex numbers in input file? */
2015 int vertex; /* current vertex being read */
2016 int new_vertex; /* new vertex being read */
2017 float weight; /* weight being read */
2018 float eweight; /* edge weight being read */
2019 int neighbor; /* neighbor of current vertex */
2020 int error_flag; /* error reading input? */
2021 int j; /* loop counters */
2022
2023 /* Read first line of input (= nvtxs, narcs, option). */
2024 /* The (decimal) digits of the option variable mean: 1's digit not zero => input
2025 edge weights 10's digit not zero => input vertex weights 100's digit not zero
2026 => include vertex numbers */
2027
2028 *start = NULL;
2029 *adjacency = NULL;
2030 *vweights = NULL;
2031 *eweights = NULL;
2032
2033 error_flag = 0;
2034 line_num = 0;
2035
2036 /* Read any leading comment lines */
2037 end_flag = 1;
2038 while (end_flag == 1) {
2039 *nvtxs = chaco_read_int(fin, &end_flag);
2040 ++line_num;
2041 }
2042 if (*nvtxs <= 0) {
2043 printf("ERROR in graph file `%s':", inname);
2044 printf(" Invalid number of vertices (%d).\n", *nvtxs);
2045 fclose(fin);
2046 return(1);
2047 }
2048
2049 narcs = chaco_read_int(fin, &end_flag);
2050 if (narcs < 0) {
2051 printf("ERROR in graph file `%s':", inname);
2052 printf(" Invalid number of expected edges (%d).\n", narcs);
2053 fclose(fin);
2054 return(1);
2055 }
2056
2057 /* Check if vertex or edge weights are used */
2058 if (!end_flag) {
2059 option = chaco_read_int(fin, &end_flag);
2060 }
2061 using_ewgts = option - 10 * (option / 10);
2062 option /= 10;
2063 using_vwgts = option - 10 * (option / 10);
2064 option /= 10;
2065 vtxnums = option - 10 * (option / 10);
2066
2067 /* Get weight info from Chaco option */
2068 (*nVwgts) = using_vwgts;
2069 (*nEwgts) = using_ewgts;
2070
2071 /* Read numbers of weights if they are specified separately */
2072 if (!end_flag && using_vwgts==1){
2073 j = chaco_read_int(fin, &end_flag);
2074 if (!end_flag) (*nVwgts) = j;
2075 }
2076 if (!end_flag && using_ewgts==1){
2077 j = chaco_read_int(fin, &end_flag);
2078 if (!end_flag) (*nEwgts) = j;
2079 }
2080
2081 /* Discard rest of line */
2082 while (!end_flag)
2083 j = chaco_read_int(fin, &end_flag);
2084
2085 /* Allocate space for rows and columns. */
2086 *start = (int *) malloc((unsigned) (*nvtxs + 1) * sizeof(int));
2087 if (narcs != 0)
2088 *adjacency = (int *) malloc((unsigned) (2 * narcs + 1) * sizeof(int));
2089 else
2090 *adjacency = NULL;
2091
2092 if (using_vwgts)
2093 *vweights = (float *) malloc((unsigned) (*nvtxs) * (*nVwgts) * sizeof(float));
2094 else
2095 *vweights = NULL;
2096
2097 if (using_ewgts)
2098 *eweights = (float *)
2099 malloc((unsigned) (2 * narcs + 1) * (*nEwgts) * sizeof(float));
2100 else
2101 *eweights = NULL;
2102
2103 adjptr = *adjacency;
2104 ewptr = *eweights;
2105
2106 sum_edges = 0;
2107 nedges = 0;
2108 (*start)[0] = 0;
2109 vertex = 0;
2110 vtx = 0;
2111 new_vertex = 1;
2112 while ((using_vwgts || vtxnums || narcs) && end_flag != -1) {
2113 ++line_num;
2114
2115 /* If multiple input lines per vertex, read vertex number. */
2116 if (vtxnums) {
2117 j = chaco_read_int(fin, &end_flag);
2118 if (end_flag) {
2119 if (vertex == *nvtxs)
2120 break;
2121 printf("ERROR in graph file `%s':", inname);
2122 printf(" no vertex number in line %d.\n", line_num);
2123 fclose(fin);
2124 return (1);
2125 }
2126 if (j != vertex && j != vertex + 1) {
2127 printf("ERROR in graph file `%s':", inname);
2128 printf(" out-of-order vertex number in line %d.\n", line_num);
2129 fclose(fin);
2130 return (1);
2131 }
2132 if (j != vertex) {
2133 new_vertex = 1;
2134 vertex = j;
2135 }
2136 else
2137 new_vertex = 0;
2138 }
2139 else
2140 vertex = ++vtx;
2141
2142 if (vertex > *nvtxs)
2143 break;
2144
2145 /* If vertices are weighted, read vertex weight. */
2146 if (using_vwgts && new_vertex) {
2147 for (j=0; j<(*nVwgts); j++){
2148 weight = chaco_read_val(fin, &end_flag);
2149 if (end_flag) {
2150 printf("ERROR in graph file `%s':", inname);
2151 printf(" not enough weights for vertex %d.\n", vertex);
2152 fclose(fin);
2153 return (1);
2154 }
2155 (*vweights)[(vertex-1)*(*nVwgts)+j] = weight;
2156 }
2157 }
2158
2159 nedge = 0;
2160
2161 /* Read number of adjacent vertex. */
2162 neighbor = chaco_read_int(fin, &end_flag);
2163
2164 while (!end_flag) {
2165 skip_flag = 0;
2166
2167 if (using_ewgts) { /* Read edge weight if it's being input. */
2168 for (j=0; j<(*nEwgts); j++){
2169 eweight = chaco_read_val(fin, &end_flag);
2170
2171 if (end_flag) {
2172 printf("ERROR in graph file `%s':", inname);
2173 printf(" not enough weights for edge (%d,%d).\n", vertex, neighbor);
2174 fclose(fin);
2175 return (1);
2176 }
2177
2178 else {
2179 *ewptr++ = eweight;
2180 }
2181 }
2182 }
2183
2184 /* Add edge to data structure. */
2185 if (!skip_flag) {
2186 if (++nedges > 2*narcs) {
2187 printf("ERROR in graph file `%s':", inname);
2188 printf(" at least %d adjacencies entered, but nedges = %d\n",
2189 nedges, narcs);
2190 fclose(fin);
2191 return (1);
2192 }
2193 *adjptr++ = neighbor;
2194 nedge++;
2195 }
2196
2197 /* Read number of next adjacent vertex. */
2198 neighbor = chaco_read_int(fin, &end_flag);
2199 }
2200
2201 sum_edges += nedge;
2202 (*start)[vertex] = sum_edges;
2203 }
2204
2205 /* Make sure there's nothing else in file. */
2206 flag = 0;
2207 while (!flag && end_flag != -1) {
2208 chaco_read_int(fin, &end_flag);
2209 if (!end_flag)
2210 flag = 1;
2211 }
2212
2213 (*start)[*nvtxs] = sum_edges;
2214
2215 if (vertex != 0) { /* Normal file was read. */
2216 if (narcs) {
2217 }
2218 else { /* no edges, but did have vertex weights or vertex numbers */
2219 free(*start);
2220 *start = NULL;
2221 if (*adjacency != NULL)
2222 free(*adjacency);
2223 *adjacency = NULL;
2224 if (*eweights != NULL)
2225 free(*eweights);
2226 *eweights = NULL;
2227 }
2228 }
2229
2230 else {
2231 /* Graph was empty */
2232 free(*start);
2233 if (*adjacency != NULL)
2234 free(*adjacency);
2235 if (*vweights != NULL)
2236 free(*vweights);
2237 if (*eweights != NULL)
2238 free(*eweights);
2239 *start = NULL;
2240 *adjacency = NULL;
2241 }
2242
2243 fclose(fin);
2244
2245 return (error_flag);
2246}
2247
2248
2249int UserInputForTests::chaco_input_geom(
2250 FILE *fingeom, /* geometry input file */
2251 const char *geomname, /* name of geometry file */
2252 int nvtxs, /* number of coordinates to read */
2253 int *igeom, /* dimensionality of geometry */
2254 double **x, /* coordinates of vertices */
2255 double **y,
2256 double **z
2257)
2258{
2259 double xc, yc, zc =0; /* first x, y, z coordinate */
2260 int nread; /* number of lines of coordinates read */
2261 int flag; /* any bad data at end of file? */
2262 int line_num; /* counts input lines in file */
2263 int end_flag; /* return conditional */
2264 int ndims; /* number of values in an input line */
2265 int i=0; /* loop counter */
2266
2267 *x = *y = *z = NULL;
2268 line_num = 0;
2269 end_flag = 1;
2270 while (end_flag == 1) {
2271 xc = chaco_read_val(fingeom, &end_flag);
2272 ++line_num;
2273 }
2274
2275 if (end_flag == -1) {
2276 printf("No values found in geometry file `%s'\n", geomname);
2277 fclose(fingeom);
2278 return (1);
2279 }
2280
2281 ndims = 1;
2282 yc = chaco_read_val(fingeom, &end_flag);
2283 if (end_flag == 0) {
2284 ndims = 2;
2285 zc = chaco_read_val(fingeom, &end_flag);
2286 if (end_flag == 0) {
2287 ndims = 3;
2288 chaco_read_val(fingeom, &end_flag);
2289 if (!end_flag) {
2290 printf("Too many values on input line of geometry file `%s'\n",
2291 geomname);
2292
2293 printf(" Maximum dimensionality is 3\n");
2294 fclose(fingeom);
2295 return (1);
2296 }
2297 }
2298 }
2299
2300 *igeom = ndims;
2301
2302 *x = (double *) malloc((unsigned) nvtxs * sizeof(double));
2303 (*x)[0] = xc;
2304 if (ndims > 1) {
2305 *y = (double *) malloc((unsigned) nvtxs * sizeof(double));
2306 (*y)[0] = yc;
2307 }
2308 if (ndims > 2) {
2309 *z = (double *) malloc((unsigned) nvtxs * sizeof(double));
2310 (*z)[0] = zc;
2311 }
2312
2313 for (nread = 1; nread < nvtxs; nread++) {
2314 ++line_num;
2315 if (ndims == 1) {
2316 i = fscanf(fingeom, "%lf", &((*x)[nread]));
2317 }
2318 else if (ndims == 2) {
2319 i = fscanf(fingeom, "%lf%lf", &((*x)[nread]), &((*y)[nread]));
2320 }
2321 else if (ndims == 3) {
2322 i = fscanf(fingeom, "%lf%lf%lf", &((*x)[nread]), &((*y)[nread]),
2323 &((*z)[nread]));
2324 }
2325
2326 if (i == EOF) {
2327 printf("Too few lines of values in geometry file; nvtxs=%d, but only %d read\n",
2328 nvtxs, nread);
2329 fclose(fingeom);
2330 return (1);
2331 }
2332 else if (i != ndims) {
2333 printf("Wrong number of values in line %d of geometry file `%s'\n",
2334 line_num, geomname);
2335 fclose(fingeom);
2336 return (1);
2337 }
2338 }
2339
2340 /* Check for spurious extra stuff in file. */
2341 flag = 0;
2342 end_flag = 0;
2343 while (!flag && end_flag != -1) {
2344 chaco_read_val(fingeom, &end_flag);
2345 if (!end_flag)
2346 flag = 1;
2347 }
2348
2349 fclose(fingeom);
2350
2351 return (0);
2352}
2353
2354// Chaco input assignments from filename.assign; copied from Zoltan
2355
2356int UserInputForTests::chaco_input_assign(
2357 FILE *finassign, /* input assignment file */
2358 const char *inassignname, /* name of input assignment file */
2359 int nvtxs, /* number of vertices to output */
2360 short *assignment) /* values to be printed */
2361{
2362 int flag; /* logical conditional */
2363 int end_flag; /* return flag from read_int() */
2364 int i, j; /* loop counter */
2365
2366 /* Get the assignment vector one line at a time, checking as you go. */
2367 /* First read past any comments at top. */
2368 end_flag = 1;
2369 while (end_flag == 1) {
2370 assignment[0] = chaco_read_int(finassign, &end_flag);
2371 }
2372
2373 if (assignment[0] < 0) {
2374 printf("ERROR: Entry %d in assignment file `%s' less than zero (%d)\n",
2375 1, inassignname, assignment[0]);
2376 fclose(finassign);
2377 return (1);
2378 }
2379
2380 if (end_flag == -1) {
2381 printf("ERROR: No values found in assignment file `%s'\n", inassignname);
2382 fclose(finassign);
2383 return (1);
2384 }
2385
2386 flag = 0;
2387 int np = this->tcomm_->getSize();
2388 if (assignment[0] >= np) flag = assignment[0];
2389 srand(this->tcomm_->getRank());
2390 for (i = 1; i < nvtxs; i++) {
2391 j = fscanf(finassign, "%hd", &(assignment[i]));
2392 if (j != 1) {
2393 printf("ERROR: Too few values in assignment file `%s'.\n", inassignname);
2394 fclose(finassign);
2395 return (1);
2396 }
2397 if (assignment[i] < 0) {
2398 printf("ERROR: Entry %d in assignment file `%s' less than zero (%d)\n",
2399 i+1, inassignname, assignment[i]);
2400 fclose(finassign);
2401 return (1);
2402 }
2403 if (assignment[i] >= np) { // warn since perhaps an error -- initial part
2404 // assignment is greater than number of processors
2405 if (assignment[i] > flag)
2406 flag = assignment[i];
2407 assignment[i] = rand() % np; // randomly assign vtx to a proc in this case
2408 }
2409 }
2410 srand(rand());
2411
2412 if (flag) {
2413 printf("WARNING: Possible error in assignment file `%s'\n",
2414 inassignname);
2415 printf(" Max assignment set (%d) greater than "
2416 "max processor rank (%d)\n", flag, np-1);
2417 printf(" Some vertices given random initial assignments\n");
2418 }
2419
2420 /* Check for spurious extra stuff in file. */
2421 flag = 0;
2422 end_flag = 0;
2423 while (!flag && end_flag != -1) {
2424 chaco_read_int(finassign, &end_flag);
2425 if (!end_flag)
2426 flag = 1;
2427 }
2428 if (flag) {
2429 printf("WARNING: Possible error in assignment file `%s'\n", inassignname);
2430 printf(" Numerical data found after expected end of file\n");
2431 }
2432
2433 fclose(finassign);
2434 return (0);
2435}
2436
2437// Pamgen Reader
2438void UserInputForTests::readPamgenMeshFile(string path, string testData)
2439{
2440#ifdef HAVE_ZOLTAN2_PAMGEN
2441 int rank = this->tcomm_->getRank();
2442 if (verbose_ && tcomm_->getRank() == 0)
2443 std::cout << "UserInputForTestsBD::readPamgenFile, Read: " << testData << std::endl;
2444
2445 size_t len;
2446 std::fstream file;
2447 int dimension;
2448 if (rank == 0){
2449 // set file name
2450 std::ostringstream meshFileName;
2451 meshFileName << path << "/" << testData << ".pmgen";
2452 // open file
2453
2454 file.open(meshFileName.str(), std::ios::in);
2455
2456 if(!file.is_open()) // may be a problem with path or filename
2457 {
2458 if(verbose_ && tcomm_->getRank() == 0)
2459 {
2460 std::cout << "Unable to open pamgen mesh: ";
2461 std::cout << meshFileName.str();
2462 std::cout <<"\nPlease check file path and name." << std::endl;
2463 }
2464 len = 0; // broadcaset 0 length ->will cause exit
2465 }else{
2466 // write to character array
2467 // get size of file
2468 file.seekg (0,file.end);
2469 len = file.tellg();
2470 file.seekg (0);
2471
2472 // get dimension
2473 dimension = 2;
2474 std::string line;
2475 while(std::getline(file,line))
2476 {
2477 if( line.find("nz") != std::string::npos ||
2478 line.find("nphi") != std::string::npos)
2479 {
2480 dimension = 3;
2481 break;
2482 }
2483 }
2484
2485 file.clear();
2486 file.seekg(0, std::ios::beg);
2487 }
2488 }
2489
2490 // broadcast the file size
2491 this->tcomm_->broadcast(0,sizeof(int), (char *)&dimension);
2492 this->tcomm_->broadcast(0,sizeof(size_t),(char *)&len);
2493 this->tcomm_->barrier();
2494
2495 if(len == 0){
2496 if(verbose_ && tcomm_->getRank() == 0)
2497 std::cout << "Pamgen Mesh file size == 0, exiting UserInputForTests early." << std::endl;
2498 return;
2499 }
2500
2501 char * file_data = new char[len+1];
2502 file_data[len] = '\0'; // critical to null terminate buffer
2503 if(rank == 0){
2504 file.read(file_data,len); // if proc 0 then read file
2505 }
2506
2507 // broadcast the file to the world
2508 this->tcomm_->broadcast(0,(int)len+1,file_data);
2509 this->tcomm_->barrier();
2510
2511 // Create the PamgenMesh
2512
2513 this->pamgen_mesh = rcp(new PamgenMesh);
2514 this->havePamgenMesh = true;
2515 pamgen_mesh->createMesh(file_data,dimension,this->tcomm_);
2516
2517 // save mesh info
2518 pamgen_mesh->storeMesh();
2519 this->tcomm_->barrier();
2520
2521 // set coordinates
2522 this->setPamgenCoordinateMV();
2523
2524 // set adjacency graph
2525 this->setPamgenAdjacencyGraph();
2526
2527 this->tcomm_->barrier();
2528 if(rank == 0) file.close();
2529 delete [] file_data;
2530#else
2531 throw std::runtime_error("Pamgen requested but Trilinos "
2532 "not built with Pamgen");
2533#endif
2534}
2535
2536#ifdef HAVE_ZOLTAN2_PAMGEN
2537void UserInputForTests::setPamgenCoordinateMV()
2538{
2539 int dimension = pamgen_mesh->num_dim;
2540 // get coordinate and point info;
2541// zlno_t numLocalPoints = pamgen_mesh->num_nodes;
2542// zgno_t numGlobalPoints = pamgen_mesh->num_nodes_global;
2543 zgno_t numelements = pamgen_mesh->num_elem;
2544 zgno_t numGlobalElements = pamgen_mesh->num_elems_global;
2545 // allocate and set an array of coordinate arrays
2546 zscalar_t **elem_coords = new zscalar_t * [dimension];
2547 for(int i = 0; i < dimension; ++i){
2548 elem_coords[i] = new zscalar_t[numelements];
2549 double *tmp = &pamgen_mesh->element_coord[i*numelements];
2550 for (int j = 0; j < numelements; j++) elem_coords[i][j] = tmp[j];
2551 }
2552
2553 // make a Tpetra map
2554 RCP<Tpetra::Map<zlno_t, zgno_t, znode_t> > mp;
2555 // mp = rcp(new map_t(numGlobalElements, numelements, 0, this->tcomm_)); // constructo 1
2556
2557// Array<zgno_t>::size_type numEltsPerProc = numelements;
2558 Array<zgno_t> elementList(numelements);
2559 for (Array<zgno_t>::size_type k = 0; k < numelements; ++k) {
2560 elementList[k] = pamgen_mesh->element_order_map[k];
2561 }
2562
2563 mp = rcp (new map_t (numGlobalElements, elementList, 0, this->tcomm_)); // constructor 2
2564
2565
2566 // make an array of array views containing the coordinate data
2567 Teuchos::Array<Teuchos::ArrayView<const zscalar_t> > coordView(dimension);
2568 for (int i = 0; i < dimension; i++){
2569 if(numelements > 0){
2570 Teuchos::ArrayView<const zscalar_t> a(elem_coords[i], numelements);
2571 coordView[i] = a;
2572 }
2573 else {
2574 Teuchos::ArrayView<const zscalar_t> a;
2575 coordView[i] = a;
2576 }
2577 }
2578
2579 // set the xyz_ multivector
2580 xyz_ = RCP<tMVector_t>(new
2581 tMVector_t(mp, coordView.view(0, dimension),
2582 dimension));
2583}
2584
2585void UserInputForTests::setPamgenAdjacencyGraph()
2586{
2587// int rank = this->tcomm_->getRank();
2588// if(rank == 0) std::cout << "Making a graph from our pamgen mesh...." << std::endl;
2589
2590 // Define Types
2591// typedef zlno_t lno_t;
2592// typedef zgno_t gno_t;
2593
2594 // get info for setting up map
2595 size_t local_nodes = (size_t)this->pamgen_mesh->num_nodes;
2596 size_t local_els = (size_t)this->pamgen_mesh->num_elem;
2597
2598 size_t global_els = (size_t)this->pamgen_mesh->num_elems_global; // global rows
2599 size_t global_nodes = (size_t)this->pamgen_mesh->num_nodes_global; //global columns
2600 // make map with global elements assigned to this mesh
2601 // make range map
2602// if(rank == 0) std::cout << "Building Rowmap: " << std::endl;
2603 RCP<const map_t> rowMap = rcp(new map_t(global_els,0,this->tcomm_));
2604 RCP<const map_t> rangeMap = rowMap;
2605
2606 // make domain map
2607 RCP<const map_t> domainMap = rcp(new map_t(global_nodes,0,this->tcomm_));
2608
2609 // Get max number of nodes per element
2610 int blks = this->pamgen_mesh->num_elem_blk;
2611 int max_nodes_per_el = 0;
2612 for(int i = 0; i < blks; i++)
2613 if (this->pamgen_mesh->nodes_per_element[i] > max_nodes_per_el)
2614 max_nodes_per_el = this->pamgen_mesh->nodes_per_element[i];
2615
2616 // make the element-node adjacency matrix
2617 Teuchos::RCP<tcrsMatrix_t> C = rcp(new tcrsMatrix_t(rowMap,max_nodes_per_el));
2618
2619 Array<zgno_t> g_el_ids(local_els);
2620 for (size_t k = 0; k < local_els; ++k) {
2621 g_el_ids[k] = pamgen_mesh->global_element_numbers[k]-1;
2622 }
2623
2624 Array<zgno_t> g_node_ids(local_nodes);
2625 for (size_t k = 0; k < local_nodes; ++k) {
2626 g_node_ids[k] = pamgen_mesh->global_node_numbers[k]-1;
2627 }
2628
2629 zlno_t el_no = 0;
2630 zscalar_t one = static_cast<zscalar_t>(1);
2631
2632// if(rank == 0) std::cout << "Writing C... " << std::endl;
2633 for(int i = 0; i < blks; i++)
2634 {
2635 int el_per_block = this->pamgen_mesh->elements[i];
2636 int nodes_per_el = this->pamgen_mesh->nodes_per_element[i];
2637 int * connect = this->pamgen_mesh->elmt_node_linkage[i];
2638
2639 for(int j = 0; j < el_per_block; j++)
2640 {
2641 const zgno_t gid = static_cast<zgno_t>(g_el_ids[el_no]);
2642 for(int k = 0; k < nodes_per_el; k++)
2643 {
2644 int g_node_i = g_node_ids[connect[j*nodes_per_el+k]-1];
2645 C->insertGlobalValues(gid,
2646 Teuchos::tuple<zgno_t>(g_node_i),
2647 Teuchos::tuple<zscalar_t>(one));
2648 }
2649 el_no++;
2650 }
2651 }
2652
2653 C->fillComplete(domainMap, rangeMap);
2654
2655
2656 // Compute product C*C'
2657// if(rank == 0) std::cout << "Compute Multiplication C... " << std::endl;
2658 RCP<tcrsMatrix_t> A = rcp(new tcrsMatrix_t(rowMap,0));
2659 Tpetra::MatrixMatrix::Multiply(*C, false, *C, true, *A);
2660
2661 // remove entries not adjacent
2662 // make graph
2663// if(rank == 0) std::cout << "Writing M_... " << std::endl;
2664 this->M_ = rcp(new tcrsMatrix_t(rowMap, A->getGlobalMaxNumRowEntries()));
2665
2666// if(rank == 0) std::cout << "\nSetting graph of connectivity..." << std::endl;
2667 Teuchos::ArrayView<const zgno_t> rowMapElementList =
2668 rowMap->getLocalElementList();
2669 for (Teuchos_Ordinal ii = 0; ii < rowMapElementList.size(); ii++)
2670 {
2671 zgno_t gid = rowMapElementList[ii];
2672 size_t numEntriesInRow = A->getNumEntriesInGlobalRow (gid);
2673 typename tcrsMatrix_t::nonconst_global_inds_host_view_type rowinds("Indices", numEntriesInRow);
2674 typename tcrsMatrix_t::nonconst_values_host_view_type rowvals("Values", numEntriesInRow);
2675
2676 // modified
2677 Array<zscalar_t> mod_rowvals;
2678 Array<zgno_t> mod_rowinds;
2679 A->getGlobalRowCopy (gid, rowinds, rowvals, numEntriesInRow);
2680 for (size_t i = 0; i < numEntriesInRow; i++) {
2681// if (rowvals[i] == 2*(this->pamgen_mesh->num_dim-1))
2682// {
2683 if (rowvals[i] >= 1)
2684 {
2685 mod_rowvals.push_back(one);
2686 mod_rowinds.push_back(rowinds[i]);
2687 }
2688 }
2689 this->M_->insertGlobalValues(gid, mod_rowinds, mod_rowvals);
2690 }
2691
2692 this->M_->fillComplete();
2694 // if(rank == 0) std::cout << "Completed M... " << std::endl;
2695
2696}
2697
2698#endif
2699
2700#endif
#define TEST_FAIL_AND_THROW(comm, ok, s)
const char param_comment
USERINPUT_FILE_FORMATS
A class that generates typical user input for testing.
@ MATRIX_MARKET
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
common code used by tests
float zscalar_t
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
keep typedefs that commonly appear in many places localized
Traits of Xpetra classes, including migration method.
map_t::node_type default_znode_t
RCP< xcrsGraph_t > getUIXpetraCrsGraph()
Tpetra::Map< zlno_t, zgno_t, znode_t > map_t
static void getUIRandomData(unsigned int seed, zlno_t length, zscalar_t min, zscalar_t max, ArrayView< ArrayRCP< zscalar_t > > data)
Generate lists of random scalars.
Tpetra::Export< zlno_t, zgno_t, znode_t > export_t
RCP< tVector_t > getUITpetraVector()
Tpetra::Import< zlno_t, zgno_t, znode_t > import_t
RCP< tMVector_t > getUIEdgeWeights()
RCP< xMVector_t > getUIXpetraMultiVector(int nvec)
bool hasInputDataType(const string &input_type)
RCP< tMVector_t > getUITpetraMultiVector(int nvec)
RCP< tMVector_t > getUICoordinates()
RCP< tcrsGraph_t > getUITpetraCrsGraph()
RCP< xVector_t > getUIXpetraVector()
RCP< xcrsMatrix_t > getUIXpetraCrsMatrix()
RCP< tcrsMatrix_t > getUITpetraCrsMatrix()
UserInputForTests(string path, string testData, const RCP< const Comm< int > > &c, bool debugInfo=false, bool distributeInput=true)
Constructor that reads in a matrix/graph from disk.
RCP< tMVector_t > getUIWeights()
static const std::string fail
Tpetra::Map map_t
GeometricGen::GeometricGenerator< zscalar_t, zlno_t, zgno_t, znode_t > geometricgen_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tcrsMatrix_t
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tMVector_t
Tpetra::Vector< zscalar_t, zlno_t, zgno_t, znode_t > tVector_t
fname
Begin.
Definition validXML.py:19
static RCP< tMVector_t > coordinates
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tMVector_t