Sparse left-looking rank-revealing QR factorization. More...
#include <SparseQR.h>
Public Types | |
typedef MatrixType::Index | Index |
typedef Matrix< Index, Dynamic, 1 > | IndexVector |
typedef _MatrixType | MatrixType |
typedef _OrderingType | OrderingType |
typedef PermutationMatrix < Dynamic, Dynamic, Index > | PermutationType |
typedef SparseMatrix< Scalar, ColMajor, Index > | QRMatrixType |
typedef MatrixType::RealScalar | RealScalar |
typedef MatrixType::Scalar | Scalar |
typedef Matrix< Scalar, Dynamic, 1 > | ScalarVector |
Public Member Functions | |
template<typename Rhs , typename Dest > | |
bool | _solve (const MatrixBase< Rhs > &B, MatrixBase< Dest > &dest) const |
void | analyzePattern (const MatrixType &mat) |
Preprocessing step of a QR factorization. | |
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. | |
ComputationInfo | info () const |
Reports whether previous computation was successful. | |
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 internal::solve_retval < SparseQR, Rhs > | solve (const MatrixBase< Rhs > &B) const |
template<typename Rhs > | |
const internal::sparse_solve_retval < SparseQR, Rhs > | solve (const SparseMatrixBase< Rhs > &B) const |
SparseQR () | |
SparseQR (const MatrixType &mat) | |
Protected Member Functions | |
void | sort_matrix_Q () |
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_isInitialized |
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 |
Friends | |
struct | SparseQR_QProduct |
struct | SparseQRMatrixQReturnType |
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 65 of file SparseQR.h.
typedef MatrixType::Index Eigen::SparseQR< _MatrixType, _OrderingType >::Index |
Definition at line 72 of file SparseQR.h.
typedef Matrix<Index, Dynamic, 1> Eigen::SparseQR< _MatrixType, _OrderingType >::IndexVector |
Definition at line 74 of file SparseQR.h.
typedef _MatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::MatrixType |
Definition at line 68 of file SparseQR.h.
typedef _OrderingType Eigen::SparseQR< _MatrixType, _OrderingType >::OrderingType |
Definition at line 69 of file SparseQR.h.
typedef PermutationMatrix<Dynamic, Dynamic, Index> Eigen::SparseQR< _MatrixType, _OrderingType >::PermutationType |
Definition at line 76 of file SparseQR.h.
typedef SparseMatrix<Scalar,ColMajor,Index> Eigen::SparseQR< _MatrixType, _OrderingType >::QRMatrixType |
Definition at line 73 of file SparseQR.h.
typedef MatrixType::RealScalar Eigen::SparseQR< _MatrixType, _OrderingType >::RealScalar |
Definition at line 71 of file SparseQR.h.
typedef MatrixType::Scalar Eigen::SparseQR< _MatrixType, _OrderingType >::Scalar |
Definition at line 70 of file SparseQR.h.
typedef Matrix<Scalar, Dynamic, 1> Eigen::SparseQR< _MatrixType, _OrderingType >::ScalarVector |
Definition at line 75 of file SparseQR.h.
Eigen::SparseQR< _MatrixType, _OrderingType >::SparseQR | ( | ) | [inline] |
Definition at line 78 of file SparseQR.h.
Eigen::SparseQR< _MatrixType, _OrderingType >::SparseQR | ( | const MatrixType & | mat | ) | [inline] |
Construct a QR factorization of the matrix mat.
Definition at line 87 of file SparseQR.h.
bool Eigen::SparseQR< _MatrixType, _OrderingType >::_solve | ( | const MatrixBase< Rhs > & | B, |
MatrixBase< Dest > & | dest | ||
) | const [inline] |
Definition at line 165 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 282 of file SparseQR.h.
Index Eigen::SparseQR< _MatrixType, _OrderingType >::cols | ( | void | ) | const [inline] |
Definition at line 112 of file SparseQR.h.
const PermutationType& Eigen::SparseQR< _MatrixType, _OrderingType >::colsPermutation | ( | ) | const [inline] |
Definition at line 152 of file SparseQR.h.
void Eigen::SparseQR< _MatrixType, _OrderingType >::compute | ( | const MatrixType & | mat | ) | [inline] |
Computes the QR factorization of the sparse matrix mat.
Definition at line 98 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 323 of file SparseQR.h.
ComputationInfo Eigen::SparseQR< _MatrixType, _OrderingType >::info | ( | ) | const [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 229 of file SparseQR.h.
std::string Eigen::SparseQR< _MatrixType, _OrderingType >::lastErrorMessage | ( | ) | const [inline] |
Definition at line 161 of file SparseQR.h.
SparseQRMatrixQReturnType<SparseQR> Eigen::SparseQR< _MatrixType, _OrderingType >::matrixQ | ( | void | ) | const [inline] |
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 146 of file SparseQR.h.
const QRMatrixType& Eigen::SparseQR< _MatrixType, _OrderingType >::matrixR | ( | void | ) | const [inline] |
Definition at line 116 of file SparseQR.h.
Index Eigen::SparseQR< _MatrixType, _OrderingType >::rank | ( | ) | const [inline] |
Definition at line 122 of file SparseQR.h.
Index Eigen::SparseQR< _MatrixType, _OrderingType >::rows | ( | void | ) | const [inline] |
Definition at line 108 of file SparseQR.h.
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 196 of file SparseQR.h.
const internal::solve_retval<SparseQR, Rhs> Eigen::SparseQR< _MatrixType, _OrderingType >::solve | ( | const MatrixBase< Rhs > & | B | ) | const [inline] |
Definition at line 207 of file SparseQR.h.
const internal::sparse_solve_retval<SparseQR, Rhs> Eigen::SparseQR< _MatrixType, _OrderingType >::solve | ( | const SparseMatrixBase< Rhs > & | B | ) | const [inline] |
Definition at line 214 of file SparseQR.h.
void Eigen::SparseQR< _MatrixType, _OrderingType >::sort_matrix_Q | ( | ) | [inline, protected] |
Definition at line 236 of file SparseQR.h.
friend struct SparseQR_QProduct [friend] |
Definition at line 267 of file SparseQR.h.
friend struct SparseQRMatrixQReturnType [friend] |
Definition at line 268 of file SparseQR.h.
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_analysisIsok [protected] |
Definition at line 248 of file SparseQR.h.
IndexVector Eigen::SparseQR< _MatrixType, _OrderingType >::m_etree [protected] |
Definition at line 262 of file SparseQR.h.
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_factorizationIsok [protected] |
Definition at line 249 of file SparseQR.h.
IndexVector Eigen::SparseQR< _MatrixType, _OrderingType >::m_firstRowElt [protected] |
Definition at line 263 of file SparseQR.h.
ScalarVector Eigen::SparseQR< _MatrixType, _OrderingType >::m_hcoeffs [protected] |
Definition at line 255 of file SparseQR.h.
ComputationInfo Eigen::SparseQR< _MatrixType, _OrderingType >::m_info [mutable, protected] |
Definition at line 250 of file SparseQR.h.
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_isEtreeOk [protected] |
Definition at line 265 of file SparseQR.h.
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_isInitialized [protected] |
Definition at line 247 of file SparseQR.h.
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_isQSorted [protected] |
Definition at line 264 of file SparseQR.h.
std::string Eigen::SparseQR< _MatrixType, _OrderingType >::m_lastError [protected] |
Definition at line 251 of file SparseQR.h.
Index Eigen::SparseQR< _MatrixType, _OrderingType >::m_nonzeropivots [protected] |
Definition at line 261 of file SparseQR.h.
PermutationType Eigen::SparseQR< _MatrixType, _OrderingType >::m_outputPerm_c [protected] |
Definition at line 258 of file SparseQR.h.
PermutationType Eigen::SparseQR< _MatrixType, _OrderingType >::m_perm_c [protected] |
Definition at line 256 of file SparseQR.h.
PermutationType Eigen::SparseQR< _MatrixType, _OrderingType >::m_pivotperm [protected] |
Definition at line 257 of file SparseQR.h.
QRMatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::m_pmat [protected] |
Definition at line 252 of file SparseQR.h.
QRMatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::m_Q [protected] |
Definition at line 254 of file SparseQR.h.
QRMatrixType Eigen::SparseQR< _MatrixType, _OrderingType >::m_R [protected] |
Definition at line 253 of file SparseQR.h.
RealScalar Eigen::SparseQR< _MatrixType, _OrderingType >::m_threshold [protected] |
Definition at line 259 of file SparseQR.h.
bool Eigen::SparseQR< _MatrixType, _OrderingType >::m_useDefaultThreshold [protected] |
Definition at line 260 of file SparseQR.h.