7 #ifndef SPECTRA_SPARSE_REGULAR_INVERSE_H
8 #define SPECTRA_SPARSE_REGULAR_INVERSE_H
11 #include <Eigen/SparseCore>
12 #include <Eigen/IterativeLinearSolvers>
36 template <
typename Scalar_,
int Uplo = Eigen::Lower,
int Flags = Eigen::ColMajor,
typename StorageIndex =
int>
46 using Index = Eigen::Index;
47 using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
48 using MapConstVec = Eigen::Map<const Vector>;
49 using MapVec = Eigen::Map<Vector>;
50 using SparseMatrix = Eigen::SparseMatrix<Scalar, Flags, StorageIndex>;
51 using ConstGenericSparseMatrix =
const Eigen::Ref<const SparseMatrix>;
53 ConstGenericSparseMatrix m_mat;
55 Eigen::ConjugateGradient<SparseMatrix> m_cg;
66 template <
typename Derived>
68 m_mat(mat), m_n(mat.
rows())
71 static_cast<int>(Derived::PlainObject::IsRowMajor) ==
static_cast<int>(SparseMatrix::IsRowMajor),
72 "SparseRegularInverse: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
74 if (mat.rows() != mat.cols())
75 throw std::invalid_argument(
"SparseRegularInverse: matrix must be square");
78 m_info = (m_cg.info() == Eigen::Success) ?
86 Index
rows()
const {
return m_n; }
90 Index
cols()
const {
return m_n; }
107 MapConstVec x(x_in, m_n);
108 MapVec y(y_out, m_n);
109 y.noalias() = m_cg.solve(x);
111 m_info = (m_cg.info() == Eigen::Success) ?
115 throw std::runtime_error(
"SparseRegularInverse: CG solver does not converge");
127 MapConstVec x(x_in, m_n);
128 MapVec y(y_out, m_n);
129 y.noalias() = m_mat.template selfadjointView<Uplo>() * x;
void solve(const Scalar *x_in, Scalar *y_out) const
SparseRegularInverse(const Eigen::SparseMatrixBase< Derived > &mat)
void perform_op(const Scalar *x_in, Scalar *y_out) const
@ Successful
Computation was successful.