Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SparseQR< MatrixType, OrderingType > Class Template Reference

Detailed Description

template<typename MatrixType, typename OrderingType>
class Eigen::SparseQR< MatrixType, OrderingType >

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.

Template Parameters
_MatrixTypeThe type of the sparse matrix A, must be a column-major SparseMatrix<>
_OrderingTypeThe fill-reducing ordering method. See the OrderingMethods module for the list of built-in and external ordering methods.
Warning
The input sparse matrix A must be in compressed mode (see SparseMatrix::makeCompressed()).

Public Member Functions

void analyzePattern (const MatrixType &mat)
 Preprocessing step of a QR factorization. More...
 
Index cols () const
 
const PermutationTypecolsPermutation () 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 QRMatrixTypematrixR () 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)
 

Constructor & Destructor Documentation

SparseQR ( const MatrixType &  mat)
inline

Construct a QR factorization of the matrix mat.

Warning
The matrix mat must be in compressed mode (see SparseMatrix::makeCompressed()).
See Also
compute()

References SparseQR< MatrixType, OrderingType >::compute().

Member Function Documentation

void analyzePattern ( const MatrixType &  mat)

Preprocessing step of a QR factorization.

Warning
The matrix mat must be in compressed mode (see SparseMatrix::makeCompressed()).

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.

Note
In this step it is assumed that there is no empty row in the matrix mat.

References Eigen::internal::coletree().

Referenced by SparseQR< MatrixType, OrderingType >::compute().

Index cols ( void  ) const
inline
Returns
the number of columns of the represented matrix.

References SparseMatrix< _Scalar, _Options, _Index >::cols().

const PermutationType& colsPermutation ( ) const
inline
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.
void compute ( const MatrixType &  mat)
inline

Computes the QR factorization of the sparse matrix mat.

Warning
The matrix mat must be in compressed mode (see SparseMatrix::makeCompressed()).
See Also
analyzePattern(), factorize()

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.

Parameters
matThe 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().

ComputationInfo info ( ) const
inline

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
See Also
iparm()
std::string lastErrorMessage ( ) const
inline
Returns
A string describing the type of error. This method is provided to ease debugging, not to handle errors.
SparseQRMatrixQReturnType<SparseQR> matrixQ ( void  ) const
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

* 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.

const QRMatrixType& matrixR ( ) const
inline
Returns
a const reference to the sparse upper triangular matrix R of the QR factorization.
Index rank ( ) const
inline
Returns
the number of non linearly dependent columns as determined by the pivoting threshold.
See Also
setPivotThreshold()
Index rows ( void  ) const
inline
Returns
the number of rows of the represented matrix.

References SparseMatrix< _Scalar, _Options, _Index >::rows().

Referenced by SparseQR< MatrixType, OrderingType >::solve().

void setPivotThreshold ( const RealScalar &  threshold)
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.

const internal::solve_retval<SparseQR, Rhs> solve ( const MatrixBase< Rhs > &  B) const
inline
Returns
the solution X of $ A X = B $ using the current decomposition of A.
See Also
compute()

References SparseQR< MatrixType, OrderingType >::rows().


The documentation for this class was generated from the following file: