Spectra  1.0.1
Header-only C++ Library for Large Scale Eigenvalue Problems
DenseSymMatProd.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_SYM_MAT_PROD_H
8 #define SPECTRA_DENSE_SYM_MAT_PROD_H
9 
10 #include <Eigen/Core>
11 
12 namespace Spectra {
13 
28 template <typename Scalar_, int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor>
30 {
31 public:
35  using Scalar = Scalar_;
36 
37 private:
38  using Index = Eigen::Index;
39  using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flags>;
40  using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
41  using MapConstVec = Eigen::Map<const Vector>;
42  using MapVec = Eigen::Map<Vector>;
43  using ConstGenericMatrix = const Eigen::Ref<const Matrix>;
44 
45  ConstGenericMatrix m_mat;
46 
47 public:
56  template <typename Derived>
57  DenseSymMatProd(const Eigen::MatrixBase<Derived>& mat) :
58  m_mat(mat)
59  {
60  static_assert(
61  static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(Matrix::IsRowMajor),
62  "DenseSymMatProd: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
63  }
64 
68  Index rows() const { return m_mat.rows(); }
72  Index cols() const { return m_mat.cols(); }
73 
80  // y_out = A * x_in
81  void perform_op(const Scalar* x_in, Scalar* y_out) const
82  {
83  MapConstVec x(x_in, m_mat.cols());
84  MapVec y(y_out, m_mat.rows());
85  y.noalias() = m_mat.template selfadjointView<Uplo>() * x;
86  }
87 
91  Matrix operator*(const Eigen::Ref<const Matrix>& mat_in) const
92  {
93  return m_mat.template selfadjointView<Uplo>() * mat_in;
94  }
95 
99  Scalar operator()(Index i, Index j) const
100  {
101  return m_mat(i, j);
102  }
103 };
104 
105 } // namespace Spectra
106 
107 #endif // SPECTRA_DENSE_SYM_MAT_PROD_H
Matrix operator*(const Eigen::Ref< const Matrix > &mat_in) const
Scalar operator()(Index i, Index j) const
void perform_op(const Scalar *x_in, Scalar *y_out) const
DenseSymMatProd(const Eigen::MatrixBase< Derived > &mat)