7 #ifndef SPECTRA_DENSE_SYM_SHIFT_SOLVE_H
8 #define SPECTRA_DENSE_SYM_SHIFT_SOLVE_H
13 #include "../LinAlg/BKLDLT.h"
14 #include "../Util/CompInfo.h"
32 template <
typename Scalar_,
int Uplo = Eigen::Lower,
int Flags = Eigen::ColMajor>
42 using Index = Eigen::Index;
43 using Matrix = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flags>;
44 using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
45 using MapConstVec = Eigen::Map<const Vector>;
46 using MapVec = Eigen::Map<Vector>;
47 using ConstGenericMatrix =
const Eigen::Ref<const Matrix>;
49 ConstGenericMatrix m_mat;
51 BKLDLT<Scalar> m_solver;
62 template <
typename Derived>
64 m_mat(mat), m_n(mat.
rows())
67 static_cast<int>(Derived::PlainObject::IsRowMajor) ==
static_cast<int>(Matrix::IsRowMajor),
68 "DenseSymShiftSolve: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
70 if (m_n != mat.cols())
71 throw std::invalid_argument(
"DenseSymShiftSolve: matrix must be square");
77 Index
rows()
const {
return m_n; }
81 Index
cols()
const {
return m_n; }
88 m_solver.compute(m_mat, Uplo, sigma);
90 throw std::invalid_argument(
"DenseSymShiftSolve: 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
DenseSymShiftSolve(const Eigen::MatrixBase< Derived > &mat)
void set_shift(const Scalar &sigma)
@ Successful
Computation was successful.