7 #ifndef SPECTRA_SPARSE_GEN_REAL_SHIFT_SOLVE_H
8 #define SPECTRA_SPARSE_GEN_REAL_SHIFT_SOLVE_H
11 #include <Eigen/SparseCore>
12 #include <Eigen/SparseLU>
30 template <
typename Scalar_,
int Flags = Eigen::ColMajor,
typename StorageIndex =
int>
40 using Index = Eigen::Index;
41 using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
42 using MapConstVec = Eigen::Map<const Vector>;
43 using MapVec = Eigen::Map<Vector>;
44 using SparseMatrix = Eigen::SparseMatrix<Scalar, Flags, StorageIndex>;
45 using ConstGenericSparseMatrix =
const Eigen::Ref<const SparseMatrix>;
47 ConstGenericSparseMatrix m_mat;
49 Eigen::SparseLU<SparseMatrix> m_solver;
59 template <
typename Derived>
61 m_mat(mat), m_n(mat.
rows())
64 static_cast<int>(Derived::PlainObject::IsRowMajor) ==
static_cast<int>(SparseMatrix::IsRowMajor),
65 "SparseGenRealShiftSolve: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
67 if (mat.rows() != mat.cols())
68 throw std::invalid_argument(
"SparseGenRealShiftSolve: matrix must be square");
74 Index
rows()
const {
return m_n; }
78 Index
cols()
const {
return m_n; }
85 SparseMatrix I(m_n, m_n);
88 m_solver.compute(m_mat - sigma * I);
89 if (m_solver.info() != Eigen::Success)
90 throw std::invalid_argument(
"SparseGenRealShiftSolve: factorization failed with the given shift");
102 MapConstVec x(x_in, m_n);
103 MapVec y(y_out, m_n);
104 y.noalias() = m_solver.solve(x);
void perform_op(const Scalar *x_in, Scalar *y_out) const
SparseGenRealShiftSolve(const Eigen::SparseMatrixBase< Derived > &mat)
void set_shift(const Scalar &sigma)