$darkmode
Eigen-unsupported  5.0.1-dev
MatrixMarketIterator.h
1 
2 // This file is part of Eigen, a lightweight C++ template library
3 // for linear algebra.
4 //
5 // Copyright (C) 2012 Desire NUENTSA WAKAM <desire.nuentsa_wakam@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_BROWSE_MATRICES_H
12 #define EIGEN_BROWSE_MATRICES_H
13 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 enum { SPD = 0x100, NonSymmetric = 0x0 };
20 
41 template <typename Scalar>
43  typedef typename NumTraits<Scalar>::Real RealScalar;
44 
45  public:
48 
49  public:
50  MatrixMarketIterator(const std::string& folder)
51  : m_sym(0), m_isvalid(false), m_matIsLoaded(false), m_hasRhs(false), m_hasrefX(false), m_folder(folder) {
52  m_folder_id = opendir(folder.c_str());
53  if (m_folder_id) Getnextvalidmatrix();
54  }
55 
57  if (m_folder_id) closedir(m_folder_id);
58  }
59 
60  inline MatrixMarketIterator& operator++() {
61  m_matIsLoaded = false;
62  m_hasrefX = false;
63  m_hasRhs = false;
64  Getnextvalidmatrix();
65  return *this;
66  }
67  inline operator bool() const { return m_isvalid; }
68 
70  inline MatrixType& matrix() {
71  // Read the matrix
72  if (m_matIsLoaded) return m_mat;
73 
74  std::string matrix_file = m_folder + "/" + m_matname + ".mtx";
75  if (!loadMarket(m_mat, matrix_file)) {
76  std::cerr << "Warning loadMarket failed when loading \"" << matrix_file << "\"" << std::endl;
77  m_matIsLoaded = false;
78  return m_mat;
79  }
80  m_matIsLoaded = true;
81 
82  if (m_sym != NonSymmetric) {
83  // Check whether we need to restore a full matrix:
84  RealScalar diag_norm = m_mat.diagonal().norm();
85  RealScalar lower_norm = m_mat.template triangularView<Lower>().norm();
86  RealScalar upper_norm = m_mat.template triangularView<Upper>().norm();
87  if (lower_norm > diag_norm && upper_norm == diag_norm) {
88  // only the lower part is stored
89  MatrixType tmp(m_mat);
90  m_mat = tmp.template selfadjointView<Lower>();
91  } else if (upper_norm > diag_norm && lower_norm == diag_norm) {
92  // only the upper part is stored
93  MatrixType tmp(m_mat);
94  m_mat = tmp.template selfadjointView<Upper>();
95  }
96  }
97  return m_mat;
98  }
99 
103  inline VectorType& rhs() {
104  // Get the right hand side
105  if (m_hasRhs) return m_rhs;
106 
107  std::string rhs_file;
108  rhs_file = m_folder + "/" + m_matname + "_b.mtx"; // The pattern is matname_b.mtx
109  m_hasRhs = Fileexists(rhs_file);
110  if (m_hasRhs) {
111  m_rhs.resize(m_mat.cols());
112  m_hasRhs = loadMarketVector(m_rhs, rhs_file);
113  }
114  if (!m_hasRhs) {
115  // Generate a random right hand side
116  if (!m_matIsLoaded) this->matrix();
117  m_refX.resize(m_mat.cols());
118  m_refX.setRandom();
119  m_rhs = m_mat * m_refX;
120  m_hasrefX = true;
121  m_hasRhs = true;
122  }
123  return m_rhs;
124  }
125 
132  inline VectorType& refX() {
133  // Check if a reference solution is provided
134  if (m_hasrefX) return m_refX;
135 
136  std::string lhs_file;
137  lhs_file = m_folder + "/" + m_matname + "_x.mtx";
138  m_hasrefX = Fileexists(lhs_file);
139  if (m_hasrefX) {
140  m_refX.resize(m_mat.cols());
141  m_hasrefX = loadMarketVector(m_refX, lhs_file);
142  } else
143  m_refX.resize(0);
144  return m_refX;
145  }
146 
147  inline std::string& matname() { return m_matname; }
148 
149  inline int sym() { return m_sym; }
150 
151  bool hasRhs() { return m_hasRhs; }
152  bool hasrefX() { return m_hasrefX; }
153  bool isFolderValid() { return bool(m_folder_id); }
154 
155  protected:
156  inline bool Fileexists(std::string file) {
157  std::ifstream file_id(file.c_str());
158  if (!file_id.good()) {
159  return false;
160  } else {
161  file_id.close();
162  return true;
163  }
164  }
165 
166  void Getnextvalidmatrix() {
167  m_isvalid = false;
168  // Here, we return with the next valid matrix in the folder
169  while ((m_curs_id = readdir(m_folder_id)) != NULL) {
170  m_isvalid = false;
171  std::string curfile;
172  curfile = m_folder + "/" + m_curs_id->d_name;
173  // Discard if it is a folder
174  if (m_curs_id->d_type == DT_DIR) continue; // FIXME This may not be available on non BSD systems
175  // struct stat st_buf;
176  // stat (curfile.c_str(), &st_buf);
177  // if (S_ISDIR(st_buf.st_mode)) continue;
178 
179  // Determine from the header if it is a matrix or a right hand side
180  bool isvector, iscomplex = false;
181  if (!getMarketHeader(curfile, m_sym, iscomplex, isvector)) continue;
182  if (isvector) continue;
183  if (!iscomplex) {
184  if (internal::is_same<Scalar, std::complex<float> >::value ||
185  internal::is_same<Scalar, std::complex<double> >::value)
186  continue;
187  }
188  if (iscomplex) {
189  if (internal::is_same<Scalar, float>::value || internal::is_same<Scalar, double>::value) continue;
190  }
191 
192  // Get the matrix name
193  std::string filename = m_curs_id->d_name;
194  m_matname = filename.substr(0, filename.length() - 4);
195 
196  // Find if the matrix is SPD
197  size_t found = m_matname.find("SPD");
198  if ((found != std::string::npos) && (m_sym != NonSymmetric)) m_sym = SPD;
199 
200  m_isvalid = true;
201  break;
202  }
203  }
204  int m_sym; // Symmetry of the matrix
205  MatrixType m_mat; // Current matrix
206  VectorType m_rhs; // Current vector
207  VectorType m_refX; // The reference solution, if exists
208  std::string m_matname; // Matrix Name
209  bool m_isvalid;
210  bool m_matIsLoaded; // Determine if the matrix has already been loaded from the file
211  bool m_hasRhs; // The right hand side exists
212  bool m_hasrefX; // A reference solution is provided
213  std::string m_folder;
214  DIR* m_folder_id;
215  struct dirent* m_curs_id;
216 };
217 
218 } // end namespace Eigen
219 
220 #endif
MatrixType & matrix()
Definition: MatrixMarketIterator.h:70
constexpr void resize(Index rows, Index cols)
const ConstDiagonalReturnType diagonal() const
Derived & setRandom(Index size)
Namespace containing all symbols from the Eigen library.
bool getMarketHeader(const std::string &filename, int &sym, bool &iscomplex, bool &isdense)
Reads the header of a matrixmarket file and determines the properties of a matrix.
Definition: MarketIO.h:122
VectorType & rhs()
Definition: MatrixMarketIterator.h:103
bool loadMarket(SparseMatrixType &mat, const std::string &filename)
Loads a sparse matrix from a matrixmarket format file.
Definition: MarketIO.h:156
Iterator to browse matrices from a specified folder.
Definition: MatrixMarketIterator.h:42
VectorType & refX()
Definition: MatrixMarketIterator.h:132
Index cols() const
bool loadMarketVector(VectorType &vec, const std::string &filename)
Same functionality as loadMarketDense, deprecated.
Definition: MarketIO.h:284