Spectra  1.0.1
Header-only C++ Library for Large Scale Eigenvalue Problems
DenseGenMatProd.h
1 // Copyright (C) 2016-2022 Yixuan Qiu <yixuan.qiu@cos.name>
2 //
3 // This Source Code Form is subject to the terms of the Mozilla
4 // Public License v. 2.0. If a copy of the MPL was not distributed
5 // with this file, You can obtain one at https://mozilla.org/MPL/2.0/.
6 
7 #ifndef SPECTRA_DENSE_GEN_MAT_PROD_H
8 #define SPECTRA_DENSE_GEN_MAT_PROD_H
9 
10 #include <Eigen/Core>
11 
12 namespace Spectra {
13 
19 
33 template <typename Scalar_, int Flags = Eigen::ColMajor>
35 {
36 public:
40  using Scalar = Scalar_;
41 
42 private:
43  using Index = Eigen::Index;
44  using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flags>;
45  using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
46  using MapConstVec = Eigen::Map<const Vector>;
47  using MapVec = Eigen::Map<Vector>;
48  using ConstGenericMatrix = const Eigen::Ref<const Matrix>;
49 
50  ConstGenericMatrix m_mat;
51 
52 public:
61  template <typename Derived>
62  DenseGenMatProd(const Eigen::MatrixBase<Derived>& mat) :
63  m_mat(mat)
64  {
65  static_assert(
66  static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(Matrix::IsRowMajor),
67  "DenseGenMatProd: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
68  }
69 
73  Index rows() const { return m_mat.rows(); }
77  Index cols() const { return m_mat.cols(); }
78 
85  // y_out = A * x_in
86  void perform_op(const Scalar* x_in, Scalar* y_out) const
87  {
88  MapConstVec x(x_in, m_mat.cols());
89  MapVec y(y_out, m_mat.rows());
90  y.noalias() = m_mat * x;
91  }
92 
96  Matrix operator*(const Eigen::Ref<const Matrix>& mat_in) const
97  {
98  return m_mat * mat_in;
99  }
100 
104  Scalar operator()(Index i, Index j) const
105  {
106  return m_mat(i, j);
107  }
108 };
109 
110 } // namespace Spectra
111 
112 #endif // SPECTRA_DENSE_GEN_MAT_PROD_H
DenseGenMatProd(const Eigen::MatrixBase< Derived > &mat)
Matrix operator*(const Eigen::Ref< const Matrix > &mat_in) const
void perform_op(const Scalar *x_in, Scalar *y_out) const
Scalar operator()(Index i, Index j) const