Teuchos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Teuchos_MatrixMarket_Raw_Graph_Adder.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef __Teuchos_MatrixMarket_Raw_Graph_Adder_hpp
11 #define __Teuchos_MatrixMarket_Raw_Graph_Adder_hpp
12 
13 #include "Teuchos_ConfigDefs.hpp"
14 #include "Teuchos_ArrayRCP.hpp"
15 #include "Teuchos_CommHelpers.hpp"
19 
20 #include <algorithm>
21 #include <fstream>
22 #include <iostream>
23 #include <iterator>
24 #include <vector>
25 #include <stdexcept>
26 
27 
28 namespace Teuchos {
29  namespace MatrixMarket {
30  namespace Raw {
53  template<class Ordinal>
54  class GraphElement {
55  public:
58  rowIndex_ (Teuchos::OrdinalTraits<Ordinal>::invalid ()),
59  colIndex_ (Teuchos::OrdinalTraits<Ordinal>::invalid ())
60  {}
61 
63  GraphElement (const Ordinal i, const Ordinal j) :
64  rowIndex_ (i), colIndex_ (j) {}
65 
67  bool operator== (const GraphElement& rhs) {
68  return rowIndex_ == rhs.rowIndex_ && colIndex_ == rhs.colIndex_;
69  }
70 
72  bool operator!= (const GraphElement& rhs) {
73  return ! (*this == rhs);
74  }
75 
77  bool operator< (const GraphElement& rhs) const {
78  if (rowIndex_ < rhs.rowIndex_)
79  return true;
80  else if (rowIndex_ > rhs.rowIndex_)
81  return false;
82  else { // equal
83  return colIndex_ < rhs.colIndex_;
84  }
85  }
86 
88  Ordinal rowIndex() const { return rowIndex_; }
89 
91  Ordinal colIndex() const { return colIndex_; }
92 
93  private:
94  Ordinal rowIndex_, colIndex_;
95  };
96 
101  template<class Ordinal>
102  std::ostream&
103  operator<< (std::ostream& out, const GraphElement<Ordinal>& elt)
104  {
105  out << elt.rowIndex () << " " << elt.colIndex ();
106  return out;
107  }
108 
128  template<class Ordinal>
129  class GraphAdder {
130  public:
131  typedef Ordinal index_type;
133  typedef typename std::vector<element_type>::size_type size_type;
134 
146  expectedNumRows_(0),
147  expectedNumCols_(0),
149  seenNumRows_(0),
150  seenNumCols_(0),
151  seenNumEntries_(0),
152  tolerant_ (true),
153  debug_ (false)
154  {}
155 
173  GraphAdder (const Ordinal expectedNumRows,
174  const Ordinal expectedNumCols,
175  const Ordinal expectedNumEntries,
176  const bool tolerant=false,
177  const bool debug=false) :
178  expectedNumRows_(expectedNumRows),
179  expectedNumCols_(expectedNumCols),
180  expectedNumEntries_(expectedNumEntries),
181  seenNumRows_(0),
182  seenNumCols_(0),
183  seenNumEntries_(0),
184  tolerant_ (tolerant),
185  debug_ (debug)
186  {}
187 
208  void
209  operator() (const Ordinal i,
210  const Ordinal j,
211  const bool countAgainstTotal=true)
212  {
213  if (! tolerant_) {
214  const bool indexPairOutOfRange = i < 1 || j < 1 ||
216 
217  TEUCHOS_TEST_FOR_EXCEPTION(indexPairOutOfRange,
218  std::invalid_argument, "Graph is " << expectedNumRows_ << " x "
219  << expectedNumCols_ << ", so entry A(" << i << "," << j
220  << ") is out of range.");
221  if (countAgainstTotal) {
223  std::invalid_argument, "Cannot add entry A(" << i << "," << j
224  << ") to graph; already have expected number "
225  "of entries " << expectedNumEntries_ << ".");
226  }
227  }
228  // i and j are 1-based indices, but we store them as 0-based.
229  elts_.push_back (element_type (i-1, j-1));
230 
231  // Keep track of the rightmost column containing a matrix
232  // entry, and the bottommost row containing a matrix entry.
233  // This gives us a lower bound for the dimensions of the
234  // graph, and a check for the reported dimensions of the
235  // graph in the Matrix Market file.
236  seenNumRows_ = std::max (seenNumRows_, i);
237  seenNumCols_ = std::max (seenNumCols_, j);
238  if (countAgainstTotal) {
239  ++seenNumEntries_;
240  }
241  }
242 
259  void
260  print (std::ostream& out, const bool doMerge, const bool replace=false)
261  {
262  if (doMerge) {
264  (replace, std::logic_error, "replace = true not implemented!");
265  //merge (replace);
266  merge ();
267  } else {
268  std::sort (elts_.begin(), elts_.end());
269  }
270  // Print out the results, delimited by newlines.
271  typedef std::ostream_iterator<element_type> iter_type;
272  std::copy (elts_.begin(), elts_.end(), iter_type (out, "\n"));
273  }
274 
288  std::pair<size_type, size_type>
289  merge ()
290  {
291  typedef typename std::vector<element_type>::iterator iter_type;
292 
293  // Start with a sorted container. GraphElement objects sort in
294  // lexicographic order of their (row, column) indices, for
295  // easy conversion to CSR format. If you expect that the
296  // elements will usually be sorted in the desired order, you
297  // can check first whether they are already sorted. We have
298  // no such expectation, so we don't even bother to spend the
299  // extra O(# entries) operations to check.
300  std::sort (elts_.begin(), elts_.end());
301 
302  // Remove duplicate elements from the sequence
303  iter_type it;
304  it = std::unique(elts_.begin(), elts_.end());
305  size_type numUnique = std::distance(elts_.begin(),it);
306  const size_type numRemoved = elts_.size() - numUnique;
307  elts_.resize( std::distance(elts_.begin(),it) );
308  elts_.resize (numUnique);
309  return std::make_pair (numUnique, numRemoved);
310  }
311 
343  void
345  size_type& numRemovedElts,
348  {
349  using Teuchos::arcp;
350  using Teuchos::ArrayRCP;
351 
352  std::pair<size_type, size_type> mergeResult = merge();
353 
354  // At this point, elts_ is already in CSR order.
355  // Now we can allocate and fill the ind array.
356  ArrayRCP<Ordinal> ind = arcp<Ordinal> (elts_.size ());
357 
358  // Number of rows in the graph.
359  const Ordinal nrows = tolerant_ ? seenNumRows_ : expectedNumRows_;
360  ArrayRCP<Ordinal> ptr = arcp<Ordinal> (nrows + 1);
361 
362  // Copy over the elements, and fill in the ptr array with
363  // offsets. Note that merge() sorted the entries by row
364  // index, so we can assume the row indices are increasing in
365  // the list of entries.
366  Ordinal curRow = 0;
367  Ordinal curInd = 0;
368  typedef typename std::vector<element_type>::const_iterator iter_type;
369  ptr[0] = 0; // ptr always has at least one entry.
370  for (iter_type it = elts_.begin(); it != elts_.end(); ++it) {
371  const Ordinal i = it->rowIndex ();
372  const Ordinal j = it->colIndex ();
373 
374  TEUCHOS_TEST_FOR_EXCEPTION(i < curRow, std::logic_error, "The "
375  "current graph entry's row index " << i << " is less then what "
376  "should be the current row index lower bound " << curRow << ".");
377  for (Ordinal k = curRow+1; k <= i; ++k) {
378  ptr[k] = curInd;
379  }
380  curRow = i;
381 
383  static_cast<size_t> (curInd) >= elts_.size (),
384  std::logic_error, "The current index " << curInd << " into ind "
385  "is >= the number of matrix entries " << elts_.size ()
386  << ".");
387  ind[curInd] = j;
388  ++curInd;
389  }
390  for (Ordinal k = curRow+1; k <= nrows; ++k) {
391  ptr[k] = curInd;
392  }
393 
394  // Assign to outputs here, to ensure the strong exception
395  // guarantee (assuming that ArrayRCP's operator= doesn't
396  // throw).
397  rowptr = ptr;
398  colind = ind;
399  numUniqueElts = mergeResult.first;
400  numRemovedElts = mergeResult.second;
401  }
402 
404  const std::vector<element_type>& getEntries() const {
405  return elts_;
406  }
407 
409  void clear() {
410  seenNumRows_ = 0;
411  seenNumCols_ = 0;
412  seenNumEntries_ = 0;
413  elts_.resize (0);
414  }
415 
419  const Ordinal numRows() const { return seenNumRows_; }
420 
424  const Ordinal numCols() const { return seenNumCols_; }
425 
426  private:
429  bool tolerant_;
430  bool debug_;
431 
433  std::vector<element_type> elts_;
434  };
435  } // namespace Raw
436  } // namespace MatrixMarket
437 } // namespace Teuchos
438 
439 #endif // #ifndef __Teuchos_MatrixMarket_Raw_Graph_Adder_hpp
bool operator==(const GraphElement &rhs)
Compare row and column indices.
To be used with Checker for &quot;raw&quot; sparse matrix input.
bool operator<(const GraphElement &rhs) const
Lexicographic order first by row index, then by column index.
void clear()
Clear all the added graph entries and reset metadata.
bool operator!=(const GraphElement &rhs)
Compare row and column indices.
const std::vector< element_type > & getEntries() const
A temporary const view of the entries of the graph.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Macro for throwing an exception with breakpointing to ease debugging.
const Ordinal numRows() const
Computed number of rows.
const Ordinal numCols() const
Computed number of columns.
Teuchos header file which uses auto-configuration information to include necessary C++ headers...
GraphElement()
Default constructor: an invalid entry of the graph.
Templated Parameter List class.
Ordinal colIndex() const
Column index (zero-based) of this GraphElement.
GraphElement(const Ordinal i, const Ordinal j)
Create a sparse graph entry at (i,j).
void operator()(const Ordinal i, const Ordinal j, const bool countAgainstTotal=true)
Add an entry to the sparse graph.
This structure defines some basic traits for the ordinal field type.
void print(std::ostream &out, const bool doMerge, const bool replace=false)
Print the sparse graph data.
std::vector< element_type > elts_
The actual matrix entries, stored as an array of structs.
GraphAdder(const Ordinal expectedNumRows, const Ordinal expectedNumCols, const Ordinal expectedNumEntries, const bool tolerant=false, const bool debug=false)
Standard constructor.
Ordinal rowIndex() const
Row index (zero-based) of this GraphElement.
void mergeAndConvertToCSR(size_type &numUniqueElts, size_type &numRemovedElts, Teuchos::ArrayRCP< Ordinal > &rowptr, Teuchos::ArrayRCP< Ordinal > &colind)
Merge duplicate elements and convert to zero-based CSR.
std::pair< size_type, size_type > merge()
Merge duplicate elements.
Reference-counted smart pointer for managing arrays.