![]() |
Eigen
5.0.1-dev+60122df6
|
#include <Eigen/src/SPQRSupport/SuiteSparseQRSupport.h>
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<> |
This class follows the sparse solver concept .
Inheritance diagram for Eigen::SPQR< MatrixType_ >:Public Member Functions | |
| cholmod_common * | cholmodCommon () const |
| Index | cols () const |
| PermutationType | colsPermutation () const |
| Get the permutation that was applied to columns of A. | |
| ComputationInfo | info () const |
| Reports whether previous computation was successful. | |
| SPQRMatrixQReturnType< SPQR > | matrixQ () const |
| Get an expression of the matrix Q. | |
| const MatrixType | matrixR () const |
| Index | rank () const |
| Index | rows () const |
| void | setPivotThreshold (const RealScalar &tol) |
| Set the tolerance tol to treat columns with 2-norm < =tol as zero. | |
| void | setSPQROrdering (int ord) |
| Set the fill-reducing ordering method to be used. | |
Public Member Functions inherited from Eigen::SparseSolverBase< SPQR< MatrixType_ > > | |
| const Solve< SPQR< MatrixType_ >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| const Solve< SPQR< MatrixType_ >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| const Solve< SPQR< MatrixType_ >, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
| const Solve< SPQR< MatrixType_ >, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
| SparseSolverBase () | |
| SparseSolverBase () | |
|
inline |
|
inline |
Get the number of columns of the input matrix.
|
inline |
Reports whether previous computation was successful.
Success if computation was successful, NumericalIssue if the sparse QR can not be computed
|
inline |
|
inline |
Gets the rank of the matrix. It should be equal to matrixQR().cols if the matrix is full-rank
|
inline |
Get the number of rows of the input matrix and the Q matrix