Path Tracer
Public Types | Public Member Functions | Protected Types | Protected Attributes | Friends | List of all members
Eigen::SparseQR< _MatrixType, _OrderingType > Class Template Reference

Sparse left-looking QR factorization with numerical column pivoting. More...

#include <SparseQR.h>

Public Types

enum  { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef _MatrixType MatrixType
 
typedef _OrderingType OrderingType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndex > QRMatrixType
 
typedef Matrix< StorageIndex, Dynamic, 1 > IndexVector
 
typedef Matrix< Scalar, Dynamic, 1 > ScalarVector
 
typedef PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
 

Public Member Functions

 SparseQR (const MatrixType &mat)
 
void compute (const MatrixType &mat)
 
void analyzePattern (const MatrixType &mat)
 Preprocessing step of a QR factorization. More...
 
void factorize (const MatrixType &mat)
 Performs the numerical QR factorization of the input matrix. More...
 
Index rows () const
 
Index cols () const
 
const QRMatrixTypematrixR () const
 
Index rank () const
 
SparseQRMatrixQReturnType< SparseQRmatrixQ () const
 
const PermutationTypecolsPermutation () const
 
std::string lastErrorMessage () const
 
template<typename Rhs , typename Dest >
bool _solve_impl (const MatrixBase< Rhs > &B, MatrixBase< Dest > &dest) const
 
void setPivotThreshold (const RealScalar &threshold)
 
template<typename Rhs >
const Solve< SparseQR, Rhs > solve (const MatrixBase< Rhs > &B) const
 
template<typename Rhs >
const Solve< SparseQR, Rhs > solve (const SparseMatrixBase< Rhs > &B) const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
void _sort_matrix_Q ()
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Types

typedef SparseSolverBase< SparseQR< _MatrixType, _OrderingType > > Base
 

Protected Attributes

bool m_analysisIsok
 
bool m_factorizationIsok
 
ComputationInfo m_info
 
std::string m_lastError
 
QRMatrixType m_pmat
 
QRMatrixType m_R
 
QRMatrixType m_Q
 
ScalarVector m_hcoeffs
 
PermutationType m_perm_c
 
PermutationType m_pivotperm
 
PermutationType m_outputPerm_c
 
RealScalar m_threshold
 
bool m_useDefaultThreshold
 
Index m_nonzeropivots
 
IndexVector m_etree
 
IndexVector m_firstRowElt
 
bool m_isQSorted
 
bool m_isEtreeOk
 
bool m_isInitialized
 

Friends

template<typename , typename >
struct SparseQR_QProduct
 

Detailed Description

template<typename _MatrixType, typename _OrderingType>
class Eigen::SparseQR< _MatrixType, _OrderingType >

Sparse left-looking QR factorization with numerical column pivoting.

This class implements a left-looking QR decomposition of sparse matrices with numerical column pivoting. 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 numerical 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.

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.

\implsparsesolverconcept

The numerical pivoting strategy and default threshold are the same as in SuiteSparse QR, and detailed in the following paper: Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011. Even though it is qualified as "rank-revealing", this strategy might fail for some rank deficient problems. When this class is used to solve linear or least-square problems it is thus strongly recommended to check the accuracy of the computed solution. If it failed, it usually helps to increase the threshold with setPivotThreshold.

Warning
The input sparse matrix A must be in compressed mode (see SparseMatrix::makeCompressed()).
For complex matrices matrixQ().transpose() will actually return the adjoint matrix.

Constructor & Destructor Documentation

◆ SparseQR()

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseQR< _MatrixType, _OrderingType >::SparseQR ( const MatrixType &  mat)
inlineexplicit

Construct a QR factorization of the matrix mat.

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

Member Function Documentation

◆ analyzePattern()

template<typename MatrixType , typename OrderingType >
void Eigen::SparseQR< MatrixType, OrderingType >::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.

◆ cols()

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseQR< _MatrixType, _OrderingType >::cols ( void  ) const
inline
Returns
the number of columns of the represented matrix.

◆ colsPermutation()

template<typename _MatrixType , typename _OrderingType >
const PermutationType& Eigen::SparseQR< _MatrixType, _OrderingType >::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.

◆ compute()

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseQR< _MatrixType, _OrderingType >::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()

◆ factorize()

template<typename MatrixType , typename OrderingType >
void Eigen::SparseQR< MatrixType, OrderingType >::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

◆ info()

template<typename _MatrixType , typename _OrderingType >
ComputationInfo Eigen::SparseQR< _MatrixType, _OrderingType >::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()

◆ lastErrorMessage()

template<typename _MatrixType , typename _OrderingType >
std::string Eigen::SparseQR< _MatrixType, _OrderingType >::lastErrorMessage ( ) const
inline
Returns
A string describing the type of error. This method is provided to ease debugging, not to handle errors.

◆ matrixQ()

template<typename _MatrixType , typename _OrderingType >
SparseQRMatrixQReturnType<SparseQR> Eigen::SparseQR< _MatrixType, _OrderingType >::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.

◆ matrixR()

template<typename _MatrixType , typename _OrderingType >
const QRMatrixType& Eigen::SparseQR< _MatrixType, _OrderingType >::matrixR ( ) const
inline
Returns
a const reference to the sparse upper triangular matrix R of the QR factorization.
Warning
The entries of the returned matrix are not sorted. This means that using it in algorithms expecting sorted entries will fail. This include random coefficient accesses (SpaseMatrix::coeff()), and coefficient-wise operations. Matrix products and triangular solves are fine though.

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

◆ rank()

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseQR< _MatrixType, _OrderingType >::rank ( ) const
inline
Returns
the number of non linearly dependent columns as determined by the pivoting threshold.
See also
setPivotThreshold()

◆ rows()

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseQR< _MatrixType, _OrderingType >::rows ( void  ) const
inline
Returns
the number of rows of the represented matrix.

◆ setPivotThreshold()

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseQR< _MatrixType, _OrderingType >::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.

◆ solve()

template<typename _MatrixType , typename _OrderingType >
template<typename Rhs >
const Solve<SparseQR, Rhs> Eigen::SparseQR< _MatrixType, _OrderingType >::solve ( const MatrixBase< Rhs > &  B) const
inline
Returns
the solution X of $ A X = B $ using the current decomposition of A.
See also
compute()

The documentation for this class was generated from the following file:
Eigen::SparseQR::matrixQ
SparseQRMatrixQReturnType< SparseQR > matrixQ() const
Definition: SparseQR.h:186