11 #ifndef EIGEN_SPARSE_MARKET_IO_H
12 #define EIGEN_SPARSE_MARKET_IO_H
23 template <
typename Scalar,
typename StorageIndex>
24 inline void GetMarketLine (
const char* line, StorageIndex& i, StorageIndex& j, Scalar& value)
26 std::stringstream sline(line);
27 sline >>
i >>
j >> value;
30 template<>
inline void GetMarketLine (
const char* line,
int& i,
int& j,
float& value)
31 { std::sscanf(line,
"%d %d %g", &
i, &
j, &value); }
33 template<>
inline void GetMarketLine (
const char* line,
int& i,
int& j,
double& value)
34 { std::sscanf(line,
"%d %d %lg", &
i, &
j, &value); }
36 template<>
inline void GetMarketLine (
const char* line,
int& i,
int& j, std::complex<float>& value)
39 template<>
inline void GetMarketLine (
const char* line,
int& i,
int& j, std::complex<double>& value)
42 template <
typename Scalar,
typename StorageIndex>
43 inline void GetMarketLine (
const char* line, StorageIndex& i, StorageIndex& j, std::complex<Scalar>& value)
45 std::stringstream sline(line);
47 sline >>
i >>
j >> valR >> valI;
48 value = std::complex<Scalar>(valR,valI);
51 template <
typename RealScalar>
52 inline void GetDenseElt (
const std::string& line, RealScalar& val)
54 std::istringstream newline(line);
58 template <
typename RealScalar>
59 inline void GetDenseElt (
const std::string& line, std::complex<RealScalar>& val)
61 RealScalar valR, valI;
62 std::istringstream newline(line);
63 newline >> valR >> valI;
64 val = std::complex<RealScalar>(valR, valI);
67 template<
typename Scalar>
70 header=
"%%MatrixMarket matrix coordinate ";
71 if(internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value)
74 if(sym ==
Symmetric) header +=
" symmetric";
75 else if (sym ==
SelfAdjoint) header +=
" Hermitian";
76 else header +=
" general";
81 if(sym ==
Symmetric) header +=
" symmetric";
82 else header +=
" general";
86 template<
typename Scalar,
typename StorageIndex>
87 inline void PutMatrixElt(Scalar value, StorageIndex row, StorageIndex col, std::ofstream& out)
89 out <<
row <<
" "<<
col <<
" " << value <<
"\n";
91 template<
typename Scalar,
typename StorageIndex>
92 inline void PutMatrixElt(std::complex<Scalar> value, StorageIndex row, StorageIndex col, std::ofstream& out)
94 out <<
row <<
" " <<
col <<
" " << value.real() <<
" " << value.imag() <<
"\n";
98 template<
typename Scalar>
101 out << value <<
"\n";
103 template<
typename Scalar>
104 inline void putDenseElt(std::complex<Scalar> value, std::ofstream& out)
106 out << value.real() <<
" " << value.imag()<<
"\n";
122 inline bool getMarketHeader(
const std::string& filename,
int& sym,
bool& iscomplex,
bool& isdense)
127 std::ifstream in(filename.c_str(),std::ios::in);
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;
154 template<
typename SparseMatrixType>
155 bool loadMarket(SparseMatrixType& mat,
const std::string& filename)
157 typedef typename SparseMatrixType::Scalar Scalar;
158 typedef typename SparseMatrixType::StorageIndex StorageIndex;
159 std::ifstream input(filename.c_str(),std::ios::in);
164 input.rdbuf()->pubsetbuf(rdbuffer, 4096);
166 const int maxBuffersize = 2048;
167 char buffer[maxBuffersize];
169 bool readsizes =
false;
172 std::vector<T> elements;
174 Index M(-1), N(-1), NNZ(-1);
176 while(input.getline(buffer, maxBuffersize))
185 std::stringstream line(buffer);
186 line >>
M >> N >> NNZ;
192 elements.reserve(NNZ);
197 StorageIndex
i(-1),
j(-1);
203 if(
i>=0 &&
j>=0 &&
i<
M &&
j<N)
206 elements.push_back(
T(
i,
j,value));
210 std::cerr <<
"Invalid read: " <<
i <<
"," <<
j <<
"\n";
216 mat.setFromTriplets(elements.begin(), elements.end());
218 std::cerr << count <<
"!=" << NNZ <<
"\n";
235 template<
typename DenseType>
238 typedef typename DenseType::Scalar Scalar;
239 std::ifstream in(filename.c_str(), std::ios::in);
248 }
while (line[0] ==
'%');
249 std::istringstream newline(line);
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);
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";
262 std::cerr<<
"Input matrix can not be resized to"<<
rows<<
" x "<<
cols<<
"as given in " << filename <<
"\n";
286 std::cerr<<
"Unable to read all elements from file " << filename <<
"\n";
295 template<
typename VectorType>
311 template<
typename SparseMatrixType>
312 bool saveMarket(
const SparseMatrixType& mat,
const std::string& filename,
int sym = 0)
314 typedef typename SparseMatrixType::Scalar Scalar;
315 typedef typename SparseMatrixType::RealScalar RealScalar;
316 std::ofstream out(filename.c_str(),std::ios::out);
320 out.flags(std::ios_base::scientific);
321 out.precision(std::numeric_limits<RealScalar>::digits10 + 2);
323 internal::putMarketHeader<Scalar>(header, sym);
324 out << header << std::endl;
328 for(
typename SparseMatrixType::InnerIterator it(
mat,
j); it; ++it)
348 template<
typename DenseType>
351 typedef typename DenseType::Scalar Scalar;
352 typedef typename DenseType::RealScalar RealScalar;
353 std::ofstream out(filename.c_str(),std::ios::out);
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";
362 out <<
"%%MatrixMarket matrix array real general\n";
377 template<
typename VectorType>
RowXpr row(Index i) const
ColXpr col(Index i) const
EIGEN_CONSTEXPR Index outerSize() const
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
constexpr void resize(Index rows, Index cols)
void putMarketHeader(std::string &header, int sym)
void PutMatrixElt(Scalar value, StorageIndex row, StorageIndex col, std::ofstream &out)
void GetDenseElt(const std::string &line, RealScalar &val)
void GetMarketLine(const char *line, StorageIndex &i, StorageIndex &j, Scalar &value)
void putDenseElt(Scalar value, std::ofstream &out)
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