MatrixMarketIterator.h
Go to the documentation of this file.
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 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 enum {
19  SPD = 0x100,
20  NonSymmetric = 0x0
21 };
22 
43 template <typename Scalar>
45 {
47  public:
50 
51  public:
52  MatrixMarketIterator(const std::string &folder)
53  : m_sym(0), m_isvalid(false), m_matIsLoaded(false), m_hasRhs(false), m_hasrefX(false), m_folder(folder)
54  {
55  m_folder_id = opendir(folder.c_str());
56  if(m_folder_id)
58  }
59 
61  {
62  if (m_folder_id) closedir(m_folder_id);
63  }
64 
66  {
67  m_matIsLoaded = false;
68  m_hasrefX = false;
69  m_hasRhs = false;
71  return *this;
72  }
73  inline operator bool() const { return m_isvalid;}
74 
76  inline MatrixType& matrix()
77  {
78  // Read the matrix
79  if (m_matIsLoaded) return m_mat;
80 
81  std::string matrix_file = m_folder + "/" + m_matname + ".mtx";
82  if ( !loadMarket(m_mat, matrix_file))
83  {
84  std::cerr << "Warning loadMarket failed when loading \"" << matrix_file << "\"" << std::endl;
85  m_matIsLoaded = false;
86  return m_mat;
87  }
88  m_matIsLoaded = true;
89 
90  if (m_sym != NonSymmetric)
91  {
92  // Check whether we need to restore a full matrix:
93  RealScalar diag_norm = m_mat.diagonal().norm();
94  RealScalar lower_norm = m_mat.template triangularView<Lower>().norm();
95  RealScalar upper_norm = m_mat.template triangularView<Upper>().norm();
96  if(lower_norm>diag_norm && upper_norm==diag_norm)
97  {
98  // only the lower part is stored
99  MatrixType tmp(m_mat);
100  m_mat = tmp.template selfadjointView<Lower>();
101  }
102  else if(upper_norm>diag_norm && lower_norm==diag_norm)
103  {
104  // only the upper part is stored
105  MatrixType tmp(m_mat);
106  m_mat = tmp.template selfadjointView<Upper>();
107  }
108  }
109  return m_mat;
110  }
111 
115  inline VectorType& rhs()
116  {
117  // Get the right hand side
118  if (m_hasRhs) return m_rhs;
119 
120  std::string rhs_file;
121  rhs_file = m_folder + "/" + m_matname + "_b.mtx"; // The pattern is matname_b.mtx
122  m_hasRhs = Fileexists(rhs_file);
123  if (m_hasRhs)
124  {
125  m_rhs.resize(m_mat.cols());
126  m_hasRhs = loadMarketVector(m_rhs, rhs_file);
127  }
128  if (!m_hasRhs)
129  {
130  // Generate a random right hand side
131  if (!m_matIsLoaded) this->matrix();
132  m_refX.resize(m_mat.cols());
133  m_refX.setRandom();
134  m_rhs = m_mat * m_refX;
135  m_hasrefX = true;
136  m_hasRhs = true;
137  }
138  return m_rhs;
139  }
140 
147  inline VectorType& refX()
148  {
149  // Check if a reference solution is provided
150  if (m_hasrefX) return m_refX;
151 
152  std::string lhs_file;
153  lhs_file = m_folder + "/" + m_matname + "_x.mtx";
154  m_hasrefX = Fileexists(lhs_file);
155  if (m_hasrefX)
156  {
157  m_refX.resize(m_mat.cols());
158  m_hasrefX = loadMarketVector(m_refX, lhs_file);
159  }
160  else
161  m_refX.resize(0);
162  return m_refX;
163  }
164 
165  inline std::string& matname() { return m_matname; }
166 
167  inline int sym() { return m_sym; }
168 
169  bool hasRhs() {return m_hasRhs; }
170  bool hasrefX() {return m_hasrefX; }
171  bool isFolderValid() { return bool(m_folder_id); }
172 
173  protected:
174 
175  inline bool Fileexists(std::string file)
176  {
177  std::ifstream file_id(file.c_str());
178  if (!file_id.good() )
179  {
180  return false;
181  }
182  else
183  {
184  file_id.close();
185  return true;
186  }
187  }
188 
190  {
191  m_isvalid = false;
192  // Here, we return with the next valid matrix in the folder
193  while ( (m_curs_id = readdir(m_folder_id)) != NULL) {
194  m_isvalid = false;
195  std::string curfile;
196  curfile = m_folder + "/" + m_curs_id->d_name;
197  // Discard if it is a folder
198  if (m_curs_id->d_type == DT_DIR) continue; //FIXME This may not be available on non BSD systems
199 // struct stat st_buf;
200 // stat (curfile.c_str(), &st_buf);
201 // if (S_ISDIR(st_buf.st_mode)) continue;
202 
203  // Determine from the header if it is a matrix or a right hand side
204  bool isvector,iscomplex=false;
205  if(!getMarketHeader(curfile,m_sym,iscomplex,isvector)) continue;
206  if(isvector) continue;
207  if (!iscomplex)
208  {
209  if(internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value)
210  continue;
211  }
212  if (iscomplex)
213  {
214  if(internal::is_same<Scalar, float>::value || internal::is_same<Scalar, double>::value)
215  continue;
216  }
217 
218 
219  // Get the matrix name
220  std::string filename = m_curs_id->d_name;
221  m_matname = filename.substr(0, filename.length()-4);
222 
223  // Find if the matrix is SPD
224  size_t found = m_matname.find("SPD");
225  if( (found!=std::string::npos) && (m_sym != NonSymmetric) )
226  m_sym = SPD;
227 
228  m_isvalid = true;
229  break;
230  }
231  }
232  int m_sym; // Symmetry of the matrix
233  MatrixType m_mat; // Current matrix
234  VectorType m_rhs; // Current vector
235  VectorType m_refX; // The reference solution, if exists
236  std::string m_matname; // Matrix Name
237  bool m_isvalid;
238  bool m_matIsLoaded; // Determine if the matrix has already been loaded from the file
239  bool m_hasRhs; // The right hand side exists
240  bool m_hasrefX; // A reference solution is provided
241  std::string m_folder;
242  DIR * m_folder_id;
243  struct dirent *m_curs_id;
244 
245 };
246 
247 } // end namespace Eigen
248 
249 #endif
Iterator to browse matrices from a specified folder.
MatrixMarketIterator & operator++()
NumTraits< Scalar >::Real RealScalar
bool Fileexists(std::string file)
SparseMatrix< Scalar, ColMajor > MatrixType
Matrix< Scalar, Dynamic, 1 > VectorType
MatrixMarketIterator(const std::string &folder)
Derived & setRandom(Index rows, Index cols)
constexpr void resize(Index rows, Index cols)
RealScalar norm() const
Index cols() const
DiagonalReturnType diagonal()
bool loadMarketVector(VectorType &vec, const std::string &filename)
Same functionality as loadMarketDense, deprecated.
Definition: MarketIO.h:296
bool loadMarket(SparseMatrixType &mat, const std::string &filename)
Loads a sparse matrix from a matrixmarket format file.
Definition: MarketIO.h:155
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
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend