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

Sparse left-looking rank-revealing QR factorization. More...

#include <SparseQR.h>

Inheritance diagram for Eigen::SparseQR< _MatrixType, _OrderingType >:
Inheritance graph
[legend]

Public Types

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

Public Member Functions

template<typename Rhs , typename Dest >
bool _solve_impl (const MatrixBase< Rhs > &B, MatrixBase< Dest > &dest) const
 
void _sort_matrix_Q ()
 
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< SparseQRmatrixQ () const
 
const QRMatrixTypematrixR () const
 
Index rank () const
 
Index rows () 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
 
 SparseQR ()
 
 SparseQR (const MatrixType &mat)
 
- Public Member Functions inherited from Eigen::SparseSolverBase< SparseQR< _MatrixType, _OrderingType > >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
SparseQR< _MatrixType, _OrderingType > & derived ()
 
const SparseQR< _MatrixType, _OrderingType > & derived () const
 
const Solve< SparseQR< _MatrixType, _OrderingType >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< SparseQR< _MatrixType, _OrderingType >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 ~SparseSolverBase ()
 

Protected Types

typedef SparseSolverBase< SparseQR< _MatrixType, _OrderingType > > Base
 

Protected Attributes

bool m_analysisIsok
 
IndexVector m_etree
 
bool m_factorizationIsok
 
IndexVector m_firstRowElt
 
ScalarVector m_hcoeffs
 
ComputationInfo m_info
 
bool m_isEtreeOk
 
bool m_isQSorted
 
std::string m_lastError
 
Index m_nonzeropivots
 
PermutationType m_outputPerm_c
 
PermutationType m_perm_c
 
PermutationType m_pivotperm
 
QRMatrixType m_pmat
 
QRMatrixType m_Q
 
QRMatrixType m_R
 
RealScalar m_threshold
 
bool m_useDefaultThreshold
 
- Protected Attributes inherited from Eigen::SparseSolverBase< SparseQR< _MatrixType, _OrderingType > >
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 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()).

Definition at line 16 of file SparseQR.h.

Member Typedef Documentation

template<typename _MatrixType , typename _OrderingType >
typedef SparseSolverBase<SparseQR<_MatrixType,_OrderingType> > Eigen::SparseQR< _MatrixType, _OrderingType >::Base
protected

Definition at line 74 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef Matrix<StorageIndex, Dynamic, 1> Eigen::SparseQR< _MatrixType, _OrderingType >::IndexVector

Definition at line 84 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef _MatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::MatrixType

Definition at line 78 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef _OrderingType Eigen::SparseQR< _MatrixType, _OrderingType >::OrderingType

Definition at line 79 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> Eigen::SparseQR< _MatrixType, _OrderingType >::PermutationType

Definition at line 86 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef SparseMatrix<Scalar,ColMajor,StorageIndex> Eigen::SparseQR< _MatrixType, _OrderingType >::QRMatrixType

Definition at line 83 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::RealScalar Eigen::SparseQR< _MatrixType, _OrderingType >::RealScalar

Definition at line 81 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::Scalar Eigen::SparseQR< _MatrixType, _OrderingType >::Scalar

Definition at line 80 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef Matrix<Scalar, Dynamic, 1> Eigen::SparseQR< _MatrixType, _OrderingType >::ScalarVector

Definition at line 85 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::StorageIndex Eigen::SparseQR< _MatrixType, _OrderingType >::StorageIndex

Definition at line 82 of file SparseQR.h.

Member Enumeration Documentation

template<typename _MatrixType , typename _OrderingType >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 88 of file SparseQR.h.

Constructor & Destructor Documentation

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseQR< _MatrixType, _OrderingType >::SparseQR ( )
inline

Definition at line 94 of file SparseQR.h.

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()

Definition at line 103 of file SparseQR.h.

Member Function Documentation

template<typename _MatrixType , typename _OrderingType >
template<typename Rhs , typename Dest >
bool Eigen::SparseQR< _MatrixType, _OrderingType >::_solve_impl ( const MatrixBase< Rhs > &  B,
MatrixBase< Dest > &  dest 
) const
inline

Definition at line 192 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseQR< _MatrixType, _OrderingType >::_sort_matrix_Q ( )
inline

Definition at line 263 of file SparseQR.h.

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.

Definition at line 307 of file SparseQR.h.

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

Definition at line 128 of file SparseQR.h.

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.

Definition at line 179 of file SparseQR.h.

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()

Definition at line 114 of file SparseQR.h.

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

Definition at line 348 of file SparseQR.h.

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()

Definition at line 255 of file SparseQR.h.

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.

Definition at line 188 of file SparseQR.h.

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.

Definition at line 173 of file SparseQR.h.

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

Definition at line 143 of file SparseQR.h.

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()

Definition at line 149 of file SparseQR.h.

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

Definition at line 124 of file SparseQR.h.

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.

Definition at line 222 of file SparseQR.h.

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()

Definition at line 233 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
template<typename Rhs >
const Solve<SparseQR, Rhs> Eigen::SparseQR< _MatrixType, _OrderingType >::solve ( const SparseMatrixBase< Rhs > &  B) const
inline

Definition at line 240 of file SparseQR.h.

Friends And Related Function Documentation

template<typename _MatrixType , typename _OrderingType >
template<typename , typename >
friend struct SparseQR_QProduct
friend

Definition at line 293 of file SparseQR.h.

Member Data Documentation

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_analysisIsok
protected

Definition at line 274 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
IndexVector Eigen::SparseQR< _MatrixType, _OrderingType >::m_etree
protected

Definition at line 288 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_factorizationIsok
protected

Definition at line 275 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
IndexVector Eigen::SparseQR< _MatrixType, _OrderingType >::m_firstRowElt
protected

Definition at line 289 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
ScalarVector Eigen::SparseQR< _MatrixType, _OrderingType >::m_hcoeffs
protected

Definition at line 281 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
ComputationInfo Eigen::SparseQR< _MatrixType, _OrderingType >::m_info
mutableprotected

Definition at line 276 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_isEtreeOk
protected

Definition at line 291 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_isQSorted
protected

Definition at line 290 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
std::string Eigen::SparseQR< _MatrixType, _OrderingType >::m_lastError
protected

Definition at line 277 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseQR< _MatrixType, _OrderingType >::m_nonzeropivots
protected

Definition at line 287 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
PermutationType Eigen::SparseQR< _MatrixType, _OrderingType >::m_outputPerm_c
protected

Definition at line 284 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
PermutationType Eigen::SparseQR< _MatrixType, _OrderingType >::m_perm_c
protected

Definition at line 282 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
PermutationType Eigen::SparseQR< _MatrixType, _OrderingType >::m_pivotperm
protected

Definition at line 283 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
QRMatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::m_pmat
protected

Definition at line 278 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
QRMatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::m_Q
protected

Definition at line 280 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
QRMatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::m_R
protected

Definition at line 279 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
RealScalar Eigen::SparseQR< _MatrixType, _OrderingType >::m_threshold
protected

Definition at line 285 of file SparseQR.h.

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_useDefaultThreshold
protected

Definition at line 286 of file SparseQR.h.


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


hebiros
Author(s): Xavier Artache , Matthew Tesch
autogenerated on Thu Sep 3 2020 04:10:21