Loading [MathJax]/extensions/MathMenu.js
Spectra 1.2.0
Header-only C++ Library for Large Scale Eigenvalue Problems
DenseGenMatProd.h
1// Copyright (C) 2016-2025 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
12namespace Spectra {
13
19
34template <typename Scalar_, int Flags = Eigen::ColMajor>
36{
37public:
41 using Scalar = Scalar_;
42
43private:
44 using Index = Eigen::Index;
45 using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flags>;
46 using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
47 using MapConstVec = Eigen::Map<const Vector>;
48 using MapVec = Eigen::Map<Vector>;
49 using ConstGenericMatrix = const Eigen::Ref<const Matrix>;
50
51 ConstGenericMatrix m_mat;
52
53public:
62 template <typename Derived>
63 DenseGenMatProd(const Eigen::MatrixBase<Derived>& mat) :
64 m_mat(mat)
65 {
66 static_assert(
67 static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(Matrix::IsRowMajor),
68 "DenseGenMatProd: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
69 }
70
74 Index rows() const { return m_mat.rows(); }
78 Index cols() const { return m_mat.cols(); }
79
86 // y_out = A * x_in
87 void perform_op(const Scalar* x_in, Scalar* y_out) const
88 {
89 MapConstVec x(x_in, m_mat.cols());
90 MapVec y(y_out, m_mat.rows());
91 y.noalias() = m_mat * x;
92 }
93
97 Matrix operator*(const Eigen::Ref<const Matrix>& mat_in) const
98 {
99 return m_mat * mat_in;
100 }
101
105 Scalar operator()(Index i, Index j) const
106 {
107 return m_mat(i, j);
108 }
109};
110
111} // namespace Spectra
112
113#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