![]() |
|
Sparse QR factorization based on SuiteSparseQR library. More...
#include <SuiteSparseQRSupport.h>
Inheritance diagram for Eigen::SPQR< MatrixType_ >:Public Types | |
| enum | { ColsAtCompileTime = Dynamic , MaxColsAtCompileTime = Dynamic } |
| typedef MatrixType_::Scalar | Scalar |
| typedef MatrixType_::RealScalar | RealScalar |
| typedef SuiteSparse_long | StorageIndex |
| typedef SparseMatrix< Scalar, ColMajor, StorageIndex > | MatrixType |
| typedef Map< PermutationMatrix< Dynamic, Dynamic, StorageIndex > > | PermutationType |
Public Member Functions | |
| SPQR () | |
| SPQR (const MatrixType_ &matrix) | |
| ~SPQR () | |
| void | SPQR_free () |
| void | compute (const MatrixType_ &matrix) |
| Index | rows () const |
| Index | cols () const |
| template<typename Rhs , typename Dest > | |
| void | _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const |
| const MatrixType | matrixR () const |
| SPQRMatrixQReturnType< SPQR > | matrixQ () const |
| Get an expression of the matrix Q. More... | |
| PermutationType | colsPermutation () const |
| Get the permutation that was applied to columns of A. More... | |
| Index | rank () const |
| void | setSPQROrdering (int ord) |
| Set the fill-reducing ordering method to be used. More... | |
| void | setPivotThreshold (const RealScalar &tol) |
| Set the tolerance tol to treat columns with 2-norm < =tol as zero. More... | |
| cholmod_common * | cholmodCommon () const |
| ComputationInfo | info () const |
| Reports whether previous computation was successful. More... | |
Public Member Functions inherited from Eigen::SparseSolverBase< SPQR< MatrixType_ > > | |
| SparseSolverBase () | |
| SparseSolverBase (SparseSolverBase &&other) | |
| ~SparseSolverBase () | |
| SPQR< MatrixType_ > & | derived () |
| const SPQR< MatrixType_ > & | derived () const |
| const Solve< SPQR< MatrixType_ >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| const Solve< SPQR< MatrixType_ >, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
| void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
Protected Types | |
| typedef SparseSolverBase< SPQR< MatrixType_ > > | Base |
Protected Attributes | |
| bool | m_analysisIsOk |
| bool | m_factorizationIsOk |
| bool | m_isRUpToDate |
| ComputationInfo | m_info |
| int | m_ordering |
| int | m_allow_tol |
| RealScalar | m_tolerance |
| cholmod_sparse * | m_cR = nullptr |
| MatrixType | m_R |
| StorageIndex * | m_E = nullptr |
| cholmod_sparse * | m_H = nullptr |
| StorageIndex * | m_HPinv = nullptr |
| cholmod_dense * | m_HTau = nullptr |
| Index | m_rank |
| cholmod_common | m_cc |
| bool | m_useDefaultThreshold |
| Index | m_rows |
Protected Attributes inherited from Eigen::SparseSolverBase< SPQR< MatrixType_ > > | |
| bool | m_isInitialized |
Friends | |
| template<typename , typename > | |
| struct | SPQR_QProduct |
Sparse QR factorization based on SuiteSparseQR library.
This class is used to perform a multithreaded and multifrontal rank-revealing QR decomposition of sparse matrices. The result is then used to solve linear leasts_square systems. Clearly, a QR factorization is returned such that A*P = Q*R where :
P is the column permutation. Use colsPermutation() to get it.
Q is the orthogonal matrix represented as Householder reflectors. Use matrixQ() to get an expression and matrixQ().transpose() to get the transpose. You can then apply it to a vector.
R is the sparse triangular factor. Use matrixQR() to get it as SparseMatrix. NOTE : The Index type of R is always SuiteSparse_long. You can get it with SPQR::Index
| MatrixType_ | The type of the sparse matrix A, must be a column-major SparseMatrix<> |
\implsparsesolverconcept
|
protected |
| typedef SparseMatrix<Scalar, ColMajor, StorageIndex> Eigen::SPQR< MatrixType_ >::MatrixType |
| typedef Map<PermutationMatrix<Dynamic, Dynamic, StorageIndex> > Eigen::SPQR< MatrixType_ >::PermutationType |
| typedef MatrixType_::RealScalar Eigen::SPQR< MatrixType_ >::RealScalar |
| typedef MatrixType_::Scalar Eigen::SPQR< MatrixType_ >::Scalar |
| typedef SuiteSparse_long Eigen::SPQR< MatrixType_ >::StorageIndex |
|
inline |
References Eigen::SPQR< MatrixType_ >::m_cc.
|
inlineexplicit |
References Eigen::SPQR< MatrixType_ >::compute(), Eigen::SPQR< MatrixType_ >::m_cc, and matrix().
|
inline |
References Eigen::SPQR< MatrixType_ >::m_cc, and Eigen::SPQR< MatrixType_ >::SPQR_free().
|
inline |
References b, Eigen::SPQR< MatrixType_ >::cols(), eigen_assert, i, Eigen::SPQR< MatrixType_ >::m_E, Eigen::SPQR< MatrixType_ >::m_info, Eigen::SparseSolverBase< SPQR< MatrixType_ > >::m_isInitialized, Eigen::SPQR< MatrixType_ >::matrixQ(), Eigen::SPQR< MatrixType_ >::matrixR(), max, Eigen::SPQR< MatrixType_ >::rank(), Eigen::DenseBase< Derived >::setZero(), Eigen::Success, and y.
|
inline |
|
inline |
Get the number of columns of the input matrix.
References Eigen::SPQR< MatrixType_ >::m_cR.
Referenced by gdb.printers._MatrixEntryIterator::__next__(), Eigen::SPQR< MatrixType_ >::_solve_impl(), gdb.printers.EigenMatrixPrinter::children(), gdb.printers.EigenSparseMatrixPrinter::children(), gdb.printers.EigenMatrixPrinter::to_string(), and gdb.printers.EigenSparseMatrixPrinter::to_string().
|
inline |
Get the permutation that was applied to columns of A.
References eigen_assert, Eigen::SPQR< MatrixType_ >::m_cR, Eigen::SPQR< MatrixType_ >::m_E, and Eigen::SparseSolverBase< SPQR< MatrixType_ > >::m_isInitialized.
|
inline |
References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), Eigen::numext::is_exactly_zero(), j, Eigen::SPQR< MatrixType_ >::m_cc, Eigen::SPQR< MatrixType_ >::m_cR, Eigen::SPQR< MatrixType_ >::m_E, Eigen::SPQR< MatrixType_ >::m_H, Eigen::SPQR< MatrixType_ >::m_HPinv, Eigen::SPQR< MatrixType_ >::m_HTau, Eigen::SPQR< MatrixType_ >::m_info, Eigen::SparseSolverBase< SPQR< MatrixType_ > >::m_isInitialized, Eigen::SPQR< MatrixType_ >::m_isRUpToDate, Eigen::SPQR< MatrixType_ >::m_ordering, Eigen::SPQR< MatrixType_ >::m_rank, Eigen::SPQR< MatrixType_ >::m_rows, Eigen::SPQR< MatrixType_ >::m_tolerance, Eigen::SPQR< MatrixType_ >::m_useDefaultThreshold, matrix(), Eigen::numext::maxi(), Eigen::SparseMatrixBase< Derived >::norm(), Eigen::NumericalIssue, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows(), Eigen::SPQR< MatrixType_ >::SPQR_free(), Eigen::Success, and Eigen::viewAsCholmod().
Referenced by Eigen::SPQR< MatrixType_ >::SPQR().
|
inline |
Reports whether previous computation was successful.
Success if computation was successful, NumericalIssue if the sparse QR can not be computed References eigen_assert, Eigen::SPQR< MatrixType_ >::m_info, and Eigen::SparseSolverBase< SPQR< MatrixType_ > >::m_isInitialized.
|
inline |
Get an expression of the matrix Q.
Referenced by Eigen::SPQR< MatrixType_ >::_solve_impl().
|
inline |
References eigen_assert, Eigen::SPQR< MatrixType_ >::m_cR, Eigen::SparseSolverBase< SPQR< MatrixType_ > >::m_isInitialized, Eigen::SPQR< MatrixType_ >::m_isRUpToDate, and Eigen::SPQR< MatrixType_ >::m_R.
Referenced by Eigen::SPQR< MatrixType_ >::_solve_impl().
|
inline |
Gets the rank of the matrix. It should be equal to matrixQR().cols if the matrix is full-rank
References eigen_assert, Eigen::SPQR< MatrixType_ >::m_cc, and Eigen::SparseSolverBase< SPQR< MatrixType_ > >::m_isInitialized.
Referenced by Eigen::SPQR< MatrixType_ >::_solve_impl().
|
inline |
Get the number of rows of the input matrix and the Q matrix
References Eigen::SPQR< MatrixType_ >::m_rows.
Referenced by gdb.printers._MatrixEntryIterator::__next__(), gdb.printers.EigenMatrixPrinter::children(), gdb.printers.EigenSparseMatrixPrinter::children(), gdb.printers.EigenMatrixPrinter::to_string(), and gdb.printers.EigenSparseMatrixPrinter::to_string().
|
inline |
Set the tolerance tol to treat columns with 2-norm < =tol as zero.
References Eigen::SPQR< MatrixType_ >::m_tolerance, and Eigen::SPQR< MatrixType_ >::m_useDefaultThreshold.
|
inline |
Set the fill-reducing ordering method to be used.
References Eigen::SPQR< MatrixType_ >::m_ordering.
|
inline |
References Eigen::SPQR< MatrixType_ >::m_cc, Eigen::SPQR< MatrixType_ >::m_cR, Eigen::SPQR< MatrixType_ >::m_E, Eigen::SPQR< MatrixType_ >::m_H, Eigen::SPQR< MatrixType_ >::m_HPinv, and Eigen::SPQR< MatrixType_ >::m_HTau.
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::~SPQR().
|
protected |
|
protected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
protected |
|
mutableprotected |
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::SPQR_free().
|
mutableprotected |
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::SPQR_free().
|
mutableprotected |
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::SPQR_free().
|
mutableprotected |
|
mutableprotected |
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::matrixR().
|
protected |
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::setSPQROrdering().
|
mutableprotected |
Referenced by Eigen::SPQR< MatrixType_ >::matrixR().
|
mutableprotected |
Referenced by Eigen::SPQR< MatrixType_ >::compute().
|
protected |
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::rows().
|
protected |
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::setPivotThreshold().
|
protected |
Referenced by Eigen::SPQR< MatrixType_ >::compute(), and Eigen::SPQR< MatrixType_ >::setPivotThreshold().