File eigen_sparse_solver.h#
-
namespace sym
-
template<typename Scalar, typename EigenSolver>
class EigenSparseSolver - #include <eigen_sparse_solver.h>
A thin wrapper around Eigen’s Sparse Solver interface for use in nonlinear solver classes like sym::LevenbergMarquardtSolver.
Can be specialized with anything satisfying the SparseSolver concept.
For example, can be used like:
using LinearSolver = sym::EigenSparseSolver<double, Eigen::CholmodDecomposition<Eigen::SparseMatrix<double>>>; using NonlinearSolver = sym::LevenbergMarquardtSolver<double, LinearSolver>; using Optimizer = sym::Optimizer<double, NonlinearSolver>; Optimizer optimizer{...};
Public Types
-
using MatrixType = Eigen::SparseMatrix<Scalar>
-
using RhsType = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>
Public Functions
-
inline bool Factorize(const MatrixType &A)
Factorize A and store internally.
- Parameters:
A – a symmetric positive definite matrix.
- Returns:
true if factorization succeeded, and false if failed.
-
template<typename Rhs>
inline RhsType Solve(const Eigen::MatrixBase<Rhs> &b) const - Returns:
x for A x = b, where x and b are dense
- Pre:
this->Factorize has already been called and succeeded.
-
template<typename Rhs>
inline void SolveInPlace(Eigen::MatrixBase<Rhs> &b) const Solves in place for x in A x = b, where x and b are dense
Eigen solvers cannot actually solve in place, so this solves, then copies back into the argument.
- Pre:
this->Factorize has already been called and succeeded.
-
inline MatrixType L() const
- Returns:
the lower triangular matrix L such that P^T * L * D * L^T * P = A, where A is the last matrix to have been factorized with this->Factorize and D is a diagonal matrix with positive diagonal entries, and P is a permutation matrix.
- Pre:
this->Factorize has already been called and succeeded.
-
inline MatrixType D() const
- Returns:
the diagonal matrix D such that P^T * L * D * L^T * P = A, where A is the last matrix to have been factorized with this->Factorize, L is lower triangular with unit diagonal, and P is a permutation matrix
- Pre:
this->Factorize has already been called and succeeded.
-
inline Eigen::PermutationMatrix<Eigen::Dynamic> Permutation() const
- Returns:
the permutation matrix P such that P^T * L * D * L^T * P = A, where A is the last matrix to have been factorized with this->Factorize, L is lower triangular with unit diagonal, and D is a diagonal matrix
- Pre:
this->Factorize has already been called and succeeded.
-
inline void AnalyzeSparsityPattern(const MatrixType &A)
Defined to satisfy interface. No analysis is needed so is a no-op.
Private Members
-
EigenSolver solver_#
-
using MatrixType = Eigen::SparseMatrix<Scalar>
-
template<typename Scalar, typename EigenSolver>