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