7 #ifndef SPECTRA_DENSE_GEN_REAL_SHIFT_SOLVE_H
8 #define SPECTRA_DENSE_GEN_REAL_SHIFT_SOLVE_H
28 template <
typename Scalar_,
int Flags = Eigen::ColMajor>
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>;
45 ConstGenericMatrix m_mat;
47 Eigen::PartialPivLU<Matrix> m_solver;
58 template <
typename Derived>
60 m_mat(mat), m_n(mat.
rows())
63 static_cast<int>(Derived::PlainObject::IsRowMajor) ==
static_cast<int>(Matrix::IsRowMajor),
64 "DenseGenRealShiftSolve: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
66 if (mat.rows() != mat.cols())
67 throw std::invalid_argument(
"DenseGenRealShiftSolve: matrix must be square");
73 Index
rows()
const {
return m_n; }
77 Index
cols()
const {
return m_n; }
84 m_solver.compute(m_mat - sigma * Matrix::Identity(m_n, m_n));
96 MapConstVec x(x_in, m_n);
98 y.noalias() = m_solver.solve(x);
DenseGenRealShiftSolve(const Eigen::MatrixBase< Derived > &mat)
void set_shift(const Scalar &sigma)
void perform_op(const Scalar *x_in, Scalar *y_out) const