7 #ifndef SPECTRA_SPARSE_CHOLESKY_H
8 #define SPECTRA_SPARSE_CHOLESKY_H
11 #include <Eigen/SparseCore>
12 #include <Eigen/SparseCholesky>
15 #include "../Util/CompInfo.h"
35 template <
typename Scalar_,
int Uplo = Eigen::Lower,
int Flags = Eigen::ColMajor,
typename StorageIndex =
int>
45 using Index = Eigen::Index;
46 using Vector = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
47 using MapConstVec = Eigen::Map<const Vector>;
48 using MapVec = Eigen::Map<Vector>;
49 using SparseMatrix = Eigen::SparseMatrix<Scalar, Flags, StorageIndex>;
52 Eigen::SimplicialLLT<SparseMatrix, Uplo> m_decomp;
63 template <
typename Derived>
68 static_cast<int>(Derived::PlainObject::IsRowMajor) ==
static_cast<int>(SparseMatrix::IsRowMajor),
69 "SparseCholesky: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
71 if (mat.rows() != mat.cols())
72 throw std::invalid_argument(
"SparseCholesky: matrix must be square");
74 m_decomp.compute(mat);
75 m_info = (m_decomp.info() == Eigen::Success) ?
83 Index
rows()
const {
return m_n; }
87 Index
cols()
const {
return m_n; }
104 MapConstVec x(x_in, m_n);
105 MapVec y(y_out, m_n);
106 y.noalias() = m_decomp.permutationP() * x;
107 m_decomp.matrixL().solveInPlace(y);
119 MapConstVec x(x_in, m_n);
120 MapVec y(y_out, m_n);
121 y.noalias() = m_decomp.matrixU().solve(x);
122 y = m_decomp.permutationPinv() * y;
SparseCholesky(const Eigen::SparseMatrixBase< Derived > &mat)
void upper_triangular_solve(const Scalar *x_in, Scalar *y_out) const
void lower_triangular_solve(const Scalar *x_in, Scalar *y_out) const
@ Successful
Computation was successful.