Sparse left-looking rank-revealing QR factorization. More...
#include <SparseQR.h>
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, StorageIndex > | PermutationType |
typedef SparseMatrix< Scalar, ColMajor, StorageIndex > | QRMatrixType |
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 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 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 |
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. |
Definition at line 16 of file SparseQR.h.
|
protected |
Definition at line 74 of file SparseQR.h.
typedef Matrix<StorageIndex, Dynamic, 1> Eigen::SparseQR< _MatrixType, _OrderingType >::IndexVector |
Definition at line 84 of file SparseQR.h.
typedef _MatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::MatrixType |
Definition at line 78 of file SparseQR.h.
typedef _OrderingType Eigen::SparseQR< _MatrixType, _OrderingType >::OrderingType |
Definition at line 79 of file SparseQR.h.
typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> Eigen::SparseQR< _MatrixType, _OrderingType >::PermutationType |
Definition at line 86 of file SparseQR.h.
typedef SparseMatrix<Scalar,ColMajor,StorageIndex> Eigen::SparseQR< _MatrixType, _OrderingType >::QRMatrixType |
Definition at line 83 of file SparseQR.h.
typedef MatrixType::RealScalar Eigen::SparseQR< _MatrixType, _OrderingType >::RealScalar |
Definition at line 81 of file SparseQR.h.
typedef MatrixType::Scalar Eigen::SparseQR< _MatrixType, _OrderingType >::Scalar |
Definition at line 80 of file SparseQR.h.
typedef Matrix<Scalar, Dynamic, 1> Eigen::SparseQR< _MatrixType, _OrderingType >::ScalarVector |
Definition at line 85 of file SparseQR.h.
typedef MatrixType::StorageIndex Eigen::SparseQR< _MatrixType, _OrderingType >::StorageIndex |
Definition at line 82 of file SparseQR.h.
anonymous enum |
Enumerator | |
---|---|
ColsAtCompileTime | |
MaxColsAtCompileTime |
Definition at line 88 of file SparseQR.h.
|
inline |
Definition at line 94 of file SparseQR.h.
|
inlineexplicit |
Construct a QR factorization of the matrix mat.
Definition at line 103 of file SparseQR.h.
|
inline |
Definition at line 192 of file SparseQR.h.
|
inline |
Definition at line 263 of file SparseQR.h.
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.
Definition at line 307 of file SparseQR.h.
|
inline |
Definition at line 128 of file SparseQR.h.
|
inline |
Definition at line 179 of file SparseQR.h.
|
inline |
Computes the QR factorization of the sparse matrix mat.
Definition at line 114 of file SparseQR.h.
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.
mat | The sparse column-major matrix |
Definition at line 348 of file SparseQR.h.
|
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 invalidDefinition at line 255 of file SparseQR.h.
|
inline |
Definition at line 188 of file SparseQR.h.
|
inline |
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.
Definition at line 173 of file SparseQR.h.
|
inline |
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:
Definition at line 143 of file SparseQR.h.
|
inline |
Definition at line 149 of file SparseQR.h.
|
inline |
Definition at line 124 of file SparseQR.h.
|
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.
|
inline |
Definition at line 233 of file SparseQR.h.
|
inline |
Definition at line 240 of file SparseQR.h.
|
friend |
Definition at line 293 of file SparseQR.h.
|
protected |
Definition at line 274 of file SparseQR.h.
|
protected |
Definition at line 288 of file SparseQR.h.
|
protected |
Definition at line 275 of file SparseQR.h.
|
protected |
Definition at line 289 of file SparseQR.h.
|
protected |
Definition at line 281 of file SparseQR.h.
|
mutableprotected |
Definition at line 276 of file SparseQR.h.
|
protected |
Definition at line 291 of file SparseQR.h.
|
protected |
Definition at line 290 of file SparseQR.h.
|
protected |
Definition at line 277 of file SparseQR.h.
|
protected |
Definition at line 287 of file SparseQR.h.
|
protected |
Definition at line 284 of file SparseQR.h.
|
protected |
Definition at line 282 of file SparseQR.h.
|
protected |
Definition at line 283 of file SparseQR.h.
|
protected |
Definition at line 278 of file SparseQR.h.
|
protected |
Definition at line 280 of file SparseQR.h.
|
protected |
Definition at line 279 of file SparseQR.h.
|
protected |
Definition at line 285 of file SparseQR.h.
|
protected |
Definition at line 286 of file SparseQR.h.