OGS
EigenMatrix.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <Eigen/Sparse>
14 #include <fstream>
15 #include <string>
16 
17 #include "EigenVector.h"
20 
21 namespace MathLib
22 {
28 class EigenMatrix final
29 {
30 public:
31  using RawMatrixType = Eigen::SparseMatrix<double, Eigen::RowMajor>;
32  using IndexType = RawMatrixType::Index;
33 
34  // TODO The matrix constructor should take num_rows and num_cols as
35  // arguments
36  // that is left for a later refactoring.
37 
41  explicit EigenMatrix(IndexType const n,
42  IndexType const n_nonzero_columns = 0)
43  : mat_(n, n)
44  {
45  if (n_nonzero_columns > 0)
46  {
47  mat_.reserve(Eigen::Matrix<IndexType, Eigen::Dynamic, 1>::Constant(
48  n, n_nonzero_columns));
49  }
50  }
51 
53  IndexType getNumberOfRows() const { return mat_.rows(); }
54 
56  IndexType getNumberOfColumns() const { return mat_.cols(); }
57 
59  static constexpr IndexType getRangeBegin() { return 0; }
60 
62  IndexType getRangeEnd() const { return getNumberOfRows(); }
63 
65  void setZero()
66  {
67  auto const N = mat_.nonZeros();
68  for (auto i = decltype(N){0}; i < N; i++)
69  {
70  mat_.valuePtr()[i] = 0;
71  }
72  // don't use mat_.setZero(). it makes a matrix uncompressed
73  }
74 
77  int setValue(IndexType row, IndexType col, double val)
78  {
79  assert(row < (IndexType)getNumberOfRows() &&
80  col < (IndexType)getNumberOfColumns());
81  mat_.coeffRef(row, col) = val;
82  return 0;
83  }
84 
87  int add(IndexType row, IndexType col, double val)
88  {
89  mat_.coeffRef(row, col) += val;
90  return 0;
91  }
92 
95  template <class T_DENSE_MATRIX>
96  void add(std::vector<IndexType> const& row_pos,
97  const T_DENSE_MATRIX& sub_matrix,
98  double fkt = 1.0)
99  {
100  this->add(row_pos, row_pos, sub_matrix, fkt);
101  }
102 
105  template <class T_DENSE_MATRIX>
106  void add(RowColumnIndices<IndexType> const& indices,
107  const T_DENSE_MATRIX& sub_matrix,
108  double fkt = 1.0)
109  {
110  this->add(indices.rows, indices.columns, sub_matrix, fkt);
111  }
112 
125  template <class T_DENSE_MATRIX>
126  void add(std::vector<IndexType> const& row_pos,
127  std::vector<IndexType> const& col_pos,
128  const T_DENSE_MATRIX& sub_matrix, double fkt = 1.0);
129 
131  double get(IndexType row, IndexType col) const
132  {
133  return mat_.coeff(row, col);
134  }
135 
137  static constexpr bool isAssembled() { return true; }
138 
140  void write(const std::string& filename) const
141  {
142  std::ofstream of(filename);
143  if (of)
144  {
145  write(of);
146  }
147  }
148 
150  void write(std::ostream& os) const
151  {
152  for (int k = 0; k < mat_.outerSize(); ++k)
153  {
154  for (RawMatrixType::InnerIterator it(mat_, k); it; ++it)
155  {
156  os << it.row() << " " << it.col() << ": " << it.value() << "\n";
157  }
158  }
159  os << std::endl;
160  }
161 
163  const RawMatrixType& getRawMatrix() const { return mat_; }
164 
165 protected:
167 };
168 
169 template <class T_DENSE_MATRIX>
170 void EigenMatrix::add(std::vector<IndexType> const& row_pos,
171  std::vector<IndexType> const& col_pos,
172  const T_DENSE_MATRIX& sub_matrix, double fkt)
173 {
174  auto const n_rows = row_pos.size();
175  auto const n_cols = col_pos.size();
176  for (auto i = decltype(n_rows){0}; i < n_rows; i++)
177  {
178  auto const row = row_pos[i];
179  for (auto j = decltype(n_cols){0}; j < n_cols; j++)
180  {
181  auto const col = col_pos[j];
182  add(row, col, fkt * sub_matrix(i, j));
183  }
184  }
185 };
186 
188 template <typename SPARSITY_PATTERN>
189 struct SetMatrixSparsity<EigenMatrix, SPARSITY_PATTERN>
190 {
193  void operator()(EigenMatrix& matrix,
194  SPARSITY_PATTERN const& sparsity_pattern) const
195  {
196  static_assert(EigenMatrix::RawMatrixType::IsRowMajor,
197  "Set matrix sparsity relies on the EigenMatrix to be in "
198  "row-major storage order.");
199 
200  assert(matrix.getNumberOfRows() ==
201  static_cast<EigenMatrix::IndexType>(sparsity_pattern.size()));
202 
203  matrix.getRawMatrix().reserve(sparsity_pattern);
204  }
205 };
206 
207 } // end namespace MathLib
RawMatrixType::Index IndexType
Definition: EigenMatrix.h:32
IndexType getNumberOfColumns() const
return the number of columns
Definition: EigenMatrix.h:56
void add(std::vector< IndexType > const &row_pos, const T_DENSE_MATRIX &sub_matrix, double fkt=1.0)
Definition: EigenMatrix.h:96
static constexpr IndexType getRangeBegin()
return a start index of the active data range
Definition: EigenMatrix.h:59
void write(const std::string &filename) const
printout this matrix for debugging
Definition: EigenMatrix.h:140
void write(std::ostream &os) const
printout this matrix for debugging
Definition: EigenMatrix.h:150
void add(RowColumnIndices< IndexType > const &indices, const T_DENSE_MATRIX &sub_matrix, double fkt=1.0)
Definition: EigenMatrix.h:106
void setZero()
reset data entries to zero.
Definition: EigenMatrix.h:65
int add(IndexType row, IndexType col, double val)
Definition: EigenMatrix.h:87
RawMatrixType mat_
Definition: EigenMatrix.h:166
const RawMatrixType & getRawMatrix() const
Definition: EigenMatrix.h:163
IndexType getNumberOfRows() const
return the number of rows
Definition: EigenMatrix.h:53
int setValue(IndexType row, IndexType col, double val)
Definition: EigenMatrix.h:77
Eigen::SparseMatrix< double, Eigen::RowMajor > RawMatrixType
Definition: EigenMatrix.h:31
RawMatrixType & getRawMatrix()
Definition: EigenMatrix.h:162
double get(IndexType row, IndexType col) const
get value. This function returns zero if the element doesn't exist.
Definition: EigenMatrix.h:131
EigenMatrix(IndexType const n, IndexType const n_nonzero_columns=0)
Definition: EigenMatrix.h:41
IndexType getRangeEnd() const
return an end index of the active data range
Definition: EigenMatrix.h:62
static constexpr bool isAssembled()
return always true, i.e. the matrix is always ready for use
Definition: EigenMatrix.h:137
LineIndex const & columns
void operator()(EigenMatrix &matrix, SPARSITY_PATTERN const &sparsity_pattern) const
Definition: EigenMatrix.h:193