MarketIO.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
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_SPARSE_MARKET_IO_H
12 #define EIGEN_SPARSE_MARKET_IO_H
13 
14 #include <iostream>
15 #include <vector>
16 
17 #include "./InternalHeaderCheck.h"
18 
19 namespace Eigen {
20 
21 namespace internal
22 {
23  template <typename Scalar, typename StorageIndex>
24  inline void GetMarketLine (const char* line, StorageIndex& i, StorageIndex& j, Scalar& value)
25  {
26  std::stringstream sline(line);
27  sline >> i >> j >> value;
28  }
29 
30  template<> inline void GetMarketLine (const char* line, int& i, int& j, float& value)
31  { std::sscanf(line, "%d %d %g", &i, &j, &value); }
32 
33  template<> inline void GetMarketLine (const char* line, int& i, int& j, double& value)
34  { std::sscanf(line, "%d %d %lg", &i, &j, &value); }
35 
36  template<> inline void GetMarketLine (const char* line, int& i, int& j, std::complex<float>& value)
37  { std::sscanf(line, "%d %d %g %g", &i, &j, &numext::real_ref(value), &numext::imag_ref(value)); }
38 
39  template<> inline void GetMarketLine (const char* line, int& i, int& j, std::complex<double>& value)
40  { std::sscanf(line, "%d %d %lg %lg", &i, &j, &numext::real_ref(value), &numext::imag_ref(value)); }
41 
42  template <typename Scalar, typename StorageIndex>
43  inline void GetMarketLine (const char* line, StorageIndex& i, StorageIndex& j, std::complex<Scalar>& value)
44  {
45  std::stringstream sline(line);
46  Scalar valR, valI;
47  sline >> i >> j >> valR >> valI;
48  value = std::complex<Scalar>(valR,valI);
49  }
50 
51  template <typename RealScalar>
52  inline void GetDenseElt (const std::string& line, RealScalar& val)
53  {
54  std::istringstream newline(line);
55  newline >> val;
56  }
57 
58  template <typename RealScalar>
59  inline void GetDenseElt (const std::string& line, std::complex<RealScalar>& val)
60  {
61  RealScalar valR, valI;
62  std::istringstream newline(line);
63  newline >> valR >> valI;
64  val = std::complex<RealScalar>(valR, valI);
65  }
66 
67  template<typename Scalar>
68  inline void putMarketHeader(std::string& header,int sym)
69  {
70  header= "%%MatrixMarket matrix coordinate ";
71  if(internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value)
72  {
73  header += " complex";
74  if(sym == Symmetric) header += " symmetric";
75  else if (sym == SelfAdjoint) header += " Hermitian";
76  else header += " general";
77  }
78  else
79  {
80  header += " real";
81  if(sym == Symmetric) header += " symmetric";
82  else header += " general";
83  }
84  }
85 
86  template<typename Scalar, typename StorageIndex>
87  inline void PutMatrixElt(Scalar value, StorageIndex row, StorageIndex col, std::ofstream& out)
88  {
89  out << row << " "<< col << " " << value << "\n";
90  }
91  template<typename Scalar, typename StorageIndex>
92  inline void PutMatrixElt(std::complex<Scalar> value, StorageIndex row, StorageIndex col, std::ofstream& out)
93  {
94  out << row << " " << col << " " << value.real() << " " << value.imag() << "\n";
95  }
96 
97 
98  template<typename Scalar>
99  inline void putDenseElt(Scalar value, std::ofstream& out)
100  {
101  out << value << "\n";
102  }
103  template<typename Scalar>
104  inline void putDenseElt(std::complex<Scalar> value, std::ofstream& out)
105  {
106  out << value.real() << " " << value.imag()<< "\n";
107  }
108 
109 } // end namespace internal
110 
111 
122 inline bool getMarketHeader(const std::string& filename, int& sym, bool& iscomplex, bool& isdense)
123 {
124  sym = 0;
125  iscomplex = false;
126  isdense = false;
127  std::ifstream in(filename.c_str(),std::ios::in);
128  if(!in)
129  return false;
130 
131  std::string line;
132  // The matrix header is always the first line in the file
133  std::getline(in, line); eigen_assert(in.good());
134 
135  std::stringstream fmtline(line);
136  std::string substr[5];
137  fmtline>> substr[0] >> substr[1] >> substr[2] >> substr[3] >> substr[4];
138  if(substr[2].compare("array") == 0) isdense = true;
139  if(substr[3].compare("complex") == 0) iscomplex = true;
140  if(substr[4].compare("symmetric") == 0) sym = Symmetric;
141  else if (substr[4].compare("Hermitian") == 0) sym = SelfAdjoint;
142 
143  return true;
144 }
154 template<typename SparseMatrixType>
155 bool loadMarket(SparseMatrixType& mat, const std::string& filename)
156 {
157  typedef typename SparseMatrixType::Scalar Scalar;
158  typedef typename SparseMatrixType::StorageIndex StorageIndex;
159  std::ifstream input(filename.c_str(),std::ios::in);
160  if(!input)
161  return false;
162 
163  char rdbuffer[4096];
164  input.rdbuf()->pubsetbuf(rdbuffer, 4096);
165 
166  const int maxBuffersize = 2048;
167  char buffer[maxBuffersize];
168 
169  bool readsizes = false;
170 
172  std::vector<T> elements;
173 
174  Index M(-1), N(-1), NNZ(-1);
175  Index count = 0;
176  while(input.getline(buffer, maxBuffersize))
177  {
178  // skip comments
179  //NOTE An appropriate test should be done on the header to get the symmetry
180  if(buffer[0]=='%')
181  continue;
182 
183  if(!readsizes)
184  {
185  std::stringstream line(buffer);
186  line >> M >> N >> NNZ;
187  if(M > 0 && N > 0)
188  {
189  readsizes = true;
190  mat.resize(M,N);
191  mat.reserve(NNZ);
192  elements.reserve(NNZ);
193  }
194  }
195  else
196  {
197  StorageIndex i(-1), j(-1);
198  Scalar value;
199  internal::GetMarketLine(buffer, i, j, value);
200 
201  i--;
202  j--;
203  if(i>=0 && j>=0 && i<M && j<N)
204  {
205  ++count;
206  elements.push_back(T(i,j,value));
207  }
208  else
209  {
210  std::cerr << "Invalid read: " << i << "," << j << "\n";
211  return false;
212  }
213  }
214  }
215 
216  mat.setFromTriplets(elements.begin(), elements.end());
217  if(count!=NNZ){
218  std::cerr << count << "!=" << NNZ << "\n";
219  return false;
220  }
221  input.close();
222  return true;
223 }
224 
225 
235 template<typename DenseType>
236 bool loadMarketDense(DenseType& mat, const std::string& filename)
237 {
238  typedef typename DenseType::Scalar Scalar;
239  std::ifstream in(filename.c_str(), std::ios::in);
240  if(!in)
241  return false;
242 
243  std::string line;
244  Index rows(0), cols(0);
245  do
246  { // Skip comments
247  std::getline(in, line); eigen_assert(in.good());
248  } while (line[0] == '%');
249  std::istringstream newline(line);
250  newline >> rows >> cols;
251 
252  bool sizes_not_positive=(rows<1 || cols<1);
253  bool wrong_input_rows = (DenseType::MaxRowsAtCompileTime != Dynamic && rows > DenseType::MaxRowsAtCompileTime) ||
254  (DenseType::RowsAtCompileTime!=Dynamic && rows!=DenseType::RowsAtCompileTime);
255  bool wrong_input_cols = (DenseType::MaxColsAtCompileTime != Dynamic && cols > DenseType::MaxColsAtCompileTime) ||
256  (DenseType::ColsAtCompileTime!=Dynamic && cols!=DenseType::ColsAtCompileTime);
257 
258  if(sizes_not_positive || wrong_input_rows || wrong_input_cols){
259  if(sizes_not_positive){
260  std::cerr<< "non-positive row or column size in file" << filename << "\n";
261  }else{
262  std::cerr<< "Input matrix can not be resized to"<<rows<<" x "<<cols<< "as given in " << filename << "\n";
263  }
264  in.close();
265  return false;
266  }
267 
268  mat.resize(rows,cols);
269  Index row = 0;
270  Index col = 0;
271  Index n=0;
272  Scalar value;
273  while ( std::getline(in, line) && (row < rows) && (col < cols)){
274  internal::GetDenseElt(line, value);
275  //matrixmarket format is column major
276  mat(row,col) = value;
277  row++;
278  if(row==rows){
279  row=0;
280  col++;
281  }
282  n++;
283  }
284  in.close();
285  if (n!=mat.size()){
286  std::cerr<< "Unable to read all elements from file " << filename << "\n";
287  return false;
288  }
289  return true;
290 }
295 template<typename VectorType>
296 bool loadMarketVector(VectorType& vec, const std::string& filename)
297 {
298  return loadMarketDense(vec, filename);
299 }
300 
311 template<typename SparseMatrixType>
312 bool saveMarket(const SparseMatrixType& mat, const std::string& filename, int sym = 0)
313 {
314  typedef typename SparseMatrixType::Scalar Scalar;
315  typedef typename SparseMatrixType::RealScalar RealScalar;
316  std::ofstream out(filename.c_str(),std::ios::out);
317  if(!out)
318  return false;
319 
320  out.flags(std::ios_base::scientific);
321  out.precision(std::numeric_limits<RealScalar>::digits10 + 2);
322  std::string header;
323  internal::putMarketHeader<Scalar>(header, sym);
324  out << header << std::endl;
325  out << mat.rows() << " " << mat.cols() << " " << mat.nonZeros() << "\n";
326  int count = 0;
327  for(int j=0; j<mat.outerSize(); ++j)
328  for(typename SparseMatrixType::InnerIterator it(mat,j); it; ++it)
329  {
330  ++ count;
331  internal::PutMatrixElt(it.value(), it.row()+1, it.col()+1, out);
332  }
333  out.close();
334  return true;
335 }
336 
337 
348 template<typename DenseType>
349 bool saveMarketDense (const DenseType& mat, const std::string& filename)
350 {
351  typedef typename DenseType::Scalar Scalar;
352  typedef typename DenseType::RealScalar RealScalar;
353  std::ofstream out(filename.c_str(),std::ios::out);
354  if(!out)
355  return false;
356 
357  out.flags(std::ios_base::scientific);
358  out.precision(std::numeric_limits<RealScalar>::digits10 + 2);
359  if(internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value)
360  out << "%%MatrixMarket matrix array complex general\n";
361  else
362  out << "%%MatrixMarket matrix array real general\n";
363  out << mat.rows() << " "<< mat.cols() << "\n";
364  for (Index i=0; i < mat.cols(); i++){
365  for (Index j=0; j < mat.rows(); j++){
366  internal::putDenseElt(mat(j,i), out);
367  }
368  }
369  out.close();
370  return true;
371 }
372 
377 template<typename VectorType>
378 bool saveMarketVector (const VectorType& vec, const std::string& filename)
379 {
380  return saveMarketDense(vec, filename);
381 }
382 
383 } // end namespace Eigen
384 
385 #endif // EIGEN_SPARSE_MARKET_IO_H
int n
int i
RowXpr row(Index i) const
ColXpr col(Index i) const
Matrix4Xd M
#define eigen_assert(x)
MatrixXf mat
EIGEN_CONSTEXPR Index outerSize() const
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
bool loadMarketVector(VectorType &vec, const std::string &filename)
Same functionality as loadMarketDense, deprecated.
Definition: MarketIO.h:296
bool saveMarket(const SparseMatrixType &mat, const std::string &filename, int sym=0)
writes a sparse Matrix to a marketmarket format file
Definition: MarketIO.h:312
bool loadMarket(SparseMatrixType &mat, const std::string &filename)
Loads a sparse matrix from a matrixmarket format file.
Definition: MarketIO.h:155
bool saveMarketVector(const VectorType &vec, const std::string &filename)
Same functionality as saveMarketDense, deprecated.
Definition: MarketIO.h:378
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
bool loadMarketDense(DenseType &mat, const std::string &filename)
Loads a dense Matrix or Vector from a matrixmarket file. If a statically sized matrix has to be parse...
Definition: MarketIO.h:236
bool saveMarketDense(const DenseType &mat, const std::string &filename)
writes a dense Matrix or vector to a marketmarket format file
Definition: MarketIO.h:349
void putMarketHeader(std::string &header, int sym)
Definition: MarketIO.h:68
void PutMatrixElt(Scalar value, StorageIndex row, StorageIndex col, std::ofstream &out)
Definition: MarketIO.h:87
void GetDenseElt(const std::string &line, RealScalar &val)
Definition: MarketIO.h:52
void GetMarketLine(const char *line, StorageIndex &i, StorageIndex &j, Scalar &value)
Definition: MarketIO.h:24
void putDenseElt(Scalar value, std::ofstream &out)
Definition: MarketIO.h:99
internal::add_const_on_value_type_t< EIGEN_MATHFUNC_RETVAL(imag_ref, Scalar) > imag_ref(const Scalar &x)
internal::add_const_on_value_type_t< EIGEN_MATHFUNC_RETVAL(real_ref, Scalar) > real_ref(const Scalar &x)
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const int Dynamic
Derived::Index cols
Derived::Index rows
std::ptrdiff_t j