![]() |
Eigen
3.2.7
|
Sparse left-looking rank-revealing QR factorization.
This class implements a left-looking rank-revealing QR decomposition of sparse matrices. When a column has a norm less than a given tolerance it is implicitly permuted to the end. The QR factorization thus obtained is given by A*P = Q*R where R is upper triangular or trapezoidal.
P is the column permutation which is the product of the fill-reducing and the rank-revealing permutations. Use colsPermutation() to get it.
Q is the orthogonal matrix represented as products of 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 or trapezoidal matrix. The later occurs when A is rank-deficient. matrixR().topLeftCorner(rank(), rank()) always returns a triangular factor of full rank.
_MatrixType | The type of the sparse matrix A, must be a column-major SparseMatrix<> |
_OrderingType | The fill-reducing ordering method. See the OrderingMethods module for the list of built-in and external ordering methods. |
Public Member Functions | |
void | analyzePattern (const MatrixType &mat) |
Preprocessing step of a QR factorization. More... | |
Index | cols () const |
const PermutationType & | colsPermutation () const |
void | compute (const MatrixType &mat) |
void | factorize (const MatrixType &mat) |
Performs the numerical QR factorization of the input matrix. More... | |
ComputationInfo | info () const |
Reports whether previous computation was successful. More... | |
std::string | lastErrorMessage () const |
SparseQRMatrixQReturnType < SparseQR > | matrixQ () const |
const QRMatrixType & | matrixR () const |
Index | rank () const |
Index | rows () const |
void | setPivotThreshold (const RealScalar &threshold) |
template<typename Rhs > | |
const internal::solve_retval < SparseQR, Rhs > | solve (const MatrixBase< Rhs > &B) const |
SparseQR (const MatrixType &mat) | |
|
inline |
Construct a QR factorization of the matrix mat.
References SparseQR< MatrixType, OrderingType >::compute().
void analyzePattern | ( | const MatrixType & | mat | ) |
Preprocessing step of a QR factorization.
In this step, the fill-reducing permutation is computed and applied to the columns of A and the column elimination tree is computed as well. Only the sparsity pattern of mat is exploited.
References Eigen::internal::coletree().
Referenced by SparseQR< MatrixType, OrderingType >::compute().
|
inline |
References SparseMatrix< _Scalar, _Options, _Index >::cols().
|
inline |
|
inline |
Computes the QR factorization of the sparse matrix mat.
References SparseQR< MatrixType, OrderingType >::analyzePattern(), and SparseQR< MatrixType, OrderingType >::factorize().
Referenced by SparseQR< MatrixType, OrderingType >::SparseQR().
void factorize | ( | const MatrixType & | mat | ) |
Performs the numerical QR factorization of the input matrix.
The function SparseQR::analyzePattern(const MatrixType&) must have been called beforehand with a matrix having the same sparsity pattern than mat.
mat | The sparse column-major matrix |
References Eigen::internal::coletree(), PlainObjectBase< Derived >::data(), Eigen::InvalidInput, PlainObjectBase< Derived >::setConstant(), PlainObjectBase< Derived >::setZero(), and Eigen::Success.
Referenced by SparseQR< MatrixType, OrderingType >::compute().
|
inline |
Reports whether previous computation was successful.
Success
if computation was successful, NumericalIssue
if the QR factorization reports a numerical problem InvalidInput
if the input matrix is invalid
|
inline |
|
inline |
\returns an expression of the matrix Q as products of sparse Householder reflectors.
The common usage of this function is to apply it to a dense matrix or vector
To get a plain SparseMatrix representation of Q:
Internally, this call simply performs a sparse product between the matrix Q and a sparse identity matrix. However, due to the fact that the sparse reflectors are stored unsorted, two transpositions are needed to sort them before performing the product.
|
inline |
|
inline |
|
inline |
References SparseMatrix< _Scalar, _Options, _Index >::rows().
Referenced by SparseQR< MatrixType, OrderingType >::solve().
|
inline |
Sets the threshold that is used to determine linearly dependent columns during the factorization.
In practice, if during the factorization the norm of the column that has to be eliminated is below this threshold, then the entire column is treated as zero, and it is moved at the end.
|
inline |
References SparseQR< MatrixType, OrderingType >::rows().