template<typename _MatrixType, typename _OrderingType>
SparseQR class
Sparse left-looking rank-revealing QR factorization.
Template parameters | |
---|---|
_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. |
Contents
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().adjoint() to get the adjoint. 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.
This class follows the sparse solver concept.
Base classes
-
template<typename Derived>class SparseSolverBase
- A base class for sparse solvers.
Constructors, destructors, conversion operators
- SparseQR(const MatrixType& mat) explicit
Public functions
- void analyzePattern(const MatrixType& mat)
- Preprocessing step of a QR factorization.
- auto cols() const -> Index
- auto colsPermutation() const -> const PermutationType&
- void compute(const MatrixType& mat)
- void factorize(const MatrixType& mat)
- Performs the numerical QR factorization of the input matrix.
- auto info() const -> ComputationInfo
- Reports whether previous computation was successful.
- auto lastErrorMessage() const -> std::string
- auto matrixQ() const -> SparseQRMatrixQReturnType<SparseQR>
- auto matrixR() const -> const QRMatrixType&
- auto rank() const -> Index
- auto rows() const -> Index
- void setPivotThreshold(const RealScalar& threshold)
-
template<typename Rhs>auto solve(const MatrixBase<Rhs>& B) const -> const Solve<SparseQR, Rhs>
Function documentation
template<typename _MatrixType, typename _OrderingType>
Eigen:: SparseQR<_MatrixType, _OrderingType>:: SparseQR(const MatrixType& mat) explicit
Construct a QR factorization of the matrix mat.
template<typename _MatrixType, typename _OrderingType>
void Eigen:: SparseQR<_MatrixType, _OrderingType>:: 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.
template<typename _MatrixType, typename _OrderingType>
const PermutationType& Eigen:: SparseQR<_MatrixType, _OrderingType>:: colsPermutation() const
Returns | a const reference to the column permutation P that was applied to A such that A*P = Q*R It is the combination of the fill-in reducing permutation and numerical column pivoting. |
---|
template<typename _MatrixType, typename _OrderingType>
void Eigen:: SparseQR<_MatrixType, _OrderingType>:: compute(const MatrixType& mat)
Computes the QR factorization of the sparse matrix mat.
template<typename _MatrixType, typename _OrderingType>
void Eigen:: SparseQR<_MatrixType, _OrderingType>:: factorize(const MatrixType& mat)
Performs the numerical QR factorization of the input matrix.
Parameters | |
---|---|
mat | The sparse column-major matrix |
The function SparseQR::
template<typename _MatrixType, typename _OrderingType>
ComputationInfo Eigen:: SparseQR<_MatrixType, _OrderingType>:: info() const
Reports whether previous computation was successful.
Returns | Success if computation was successful, NumericalIssue if the QR factorization reports a numerical problem InvalidInput if the input matrix is invalid |
---|
template<typename _MatrixType, typename _OrderingType>
std::string Eigen:: SparseQR<_MatrixType, _OrderingType>:: lastErrorMessage() const
Returns | A string describing the type of error. This method is provided to ease debugging, not to handle errors. |
---|
template<typename _MatrixType, typename _OrderingType>
SparseQRMatrixQReturnType<SparseQR> Eigen:: SparseQR<_MatrixType, _OrderingType>:: matrixQ() const
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 VectorXd B1, B2; // Initialize B1 B2 = matrixQ() * B1; |
---|
To get a plain SparseMatrix representation of Q:
SparseMatrix<double> Q; Q = SparseQR<SparseMatrix<double> >(A).matrixQ();
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.
template<typename _MatrixType, typename _OrderingType>
const QRMatrixType& Eigen:: SparseQR<_MatrixType, _OrderingType>:: matrixR() const
Returns | a const reference to the sparse upper triangular matrix R of the QR factorization. |
---|
To sort the entries, you can assign it to a row-major matrix, and if a column-major matrix is required, you can copy it again:
SparseMatrix<double> R = qr.matrixR(); // column-major, not sorted! SparseMatrix<double,RowMajor> Rr = qr.matrixR(); // row-major, sorted SparseMatrix<double> Rc = Rr; // column-major, sorted
template<typename _MatrixType, typename _OrderingType>
void Eigen:: SparseQR<_MatrixType, _OrderingType>:: setPivotThreshold(const RealScalar& threshold)
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.
template<typename _MatrixType, typename _OrderingType>
template<typename Rhs>
const Solve<SparseQR, Rhs> Eigen:: SparseQR<_MatrixType, _OrderingType>:: solve(const MatrixBase<Rhs>& B) const
Returns | the solution X of using the current decomposition of A. |
---|