MOAB
4.9.3pre
|
Sparse QR factorization based on SuiteSparseQR library. More...
#include <SuiteSparseQRSupport.h>
Public Types | |
enum | { ColsAtCompileTime = Dynamic, MaxColsAtCompileTime = Dynamic } |
typedef _MatrixType::Scalar | Scalar |
typedef _MatrixType::RealScalar | RealScalar |
typedef SuiteSparse_long | StorageIndex |
typedef SparseMatrix< Scalar, ColMajor, StorageIndex > | MatrixType |
typedef Map< PermutationMatrix < Dynamic, Dynamic, StorageIndex > > | PermutationType |
Public Member Functions | |
SPQR () | |
SPQR (const _MatrixType &matrix) | |
~SPQR () | |
void | SPQR_free () |
void | compute (const _MatrixType &matrix) |
Index | rows () const |
Index | cols () const |
template<typename Rhs , typename Dest > | |
void | _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const |
const MatrixType | matrixR () const |
SPQRMatrixQReturnType< SPQR > | matrixQ () const |
Get an expression of the matrix Q. | |
PermutationType | colsPermutation () const |
Get the permutation that was applied to columns of A. | |
Index | rank () const |
void | setSPQROrdering (int ord) |
Set the fill-reducing ordering method to be used. | |
void | setPivotThreshold (const RealScalar &tol) |
Set the tolerance tol to treat columns with 2-norm < =tol as zero. | |
cholmod_common * | cholmodCommon () const |
ComputationInfo | info () const |
Reports whether previous computation was successful. | |
Protected Types | |
typedef SparseSolverBase< SPQR < _MatrixType > > | Base |
Protected Attributes | |
bool | m_analysisIsOk |
bool | m_factorizationIsOk |
bool | m_isRUpToDate |
ComputationInfo | m_info |
int | m_ordering |
int | m_allow_tol |
RealScalar | m_tolerance |
cholmod_sparse * | m_cR |
MatrixType | m_R |
StorageIndex * | m_E |
cholmod_sparse * | m_H |
StorageIndex * | m_HPinv |
cholmod_dense * | m_HTau |
Index | m_rank |
cholmod_common | m_cc |
bool | m_useDefaultThreshold |
Friends | |
struct | SPQR_QProduct |
Sparse QR factorization based on SuiteSparseQR library.
This class is used to perform a multithreaded and multifrontal rank-revealing QR decomposition of sparse matrices. The result is then used to solve linear leasts_square systems. Clearly, a QR factorization is returned such that A*P = Q*R where :
P is the column permutation. Use colsPermutation() to get it.
Q is the orthogonal matrix represented as 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 factor. Use matrixQR() to get it as SparseMatrix. NOTE : The Index type of R is always SuiteSparse_long. You can get it with SPQR::Index
_MatrixType | The type of the sparse matrix A, must be a column-major SparseMatrix<> |
Definition at line 60 of file SuiteSparseQRSupport.h.
typedef SparseSolverBase<SPQR<_MatrixType> > Eigen::SPQR< _MatrixType >::Base [protected] |
Definition at line 63 of file SuiteSparseQRSupport.h.
typedef SparseMatrix<Scalar, ColMajor, StorageIndex> Eigen::SPQR< _MatrixType >::MatrixType |
Definition at line 69 of file SuiteSparseQRSupport.h.
typedef Map<PermutationMatrix<Dynamic, Dynamic, StorageIndex> > Eigen::SPQR< _MatrixType >::PermutationType |
Definition at line 70 of file SuiteSparseQRSupport.h.
typedef _MatrixType::RealScalar Eigen::SPQR< _MatrixType >::RealScalar |
Definition at line 67 of file SuiteSparseQRSupport.h.
typedef _MatrixType::Scalar Eigen::SPQR< _MatrixType >::Scalar |
Definition at line 66 of file SuiteSparseQRSupport.h.
typedef SuiteSparse_long Eigen::SPQR< _MatrixType >::StorageIndex |
Definition at line 68 of file SuiteSparseQRSupport.h.
anonymous enum |
Definition at line 71 of file SuiteSparseQRSupport.h.
{ ColsAtCompileTime = Dynamic, MaxColsAtCompileTime = Dynamic };
Eigen::SPQR< _MatrixType >::SPQR | ( | ) | [inline] |
Definition at line 76 of file SuiteSparseQRSupport.h.
: m_ordering(SPQR_ORDERING_DEFAULT), m_allow_tol(SPQR_DEFAULT_TOL), m_tolerance (NumTraits<Scalar>::epsilon()), m_useDefaultThreshold(true) { cholmod_l_start(&m_cc); }
Eigen::SPQR< _MatrixType >::SPQR | ( | const _MatrixType & | matrix | ) | [inline, explicit] |
Definition at line 82 of file SuiteSparseQRSupport.h.
: m_ordering(SPQR_ORDERING_DEFAULT), m_allow_tol(SPQR_DEFAULT_TOL), m_tolerance (NumTraits<Scalar>::epsilon()), m_useDefaultThreshold(true) { cholmod_l_start(&m_cc); compute(matrix); }
Eigen::SPQR< _MatrixType >::~SPQR | ( | ) | [inline] |
Definition at line 89 of file SuiteSparseQRSupport.h.
void Eigen::SPQR< _MatrixType >::_solve_impl | ( | const MatrixBase< Rhs > & | b, |
MatrixBase< Dest > & | dest | ||
) | const [inline] |
Definition at line 150 of file SuiteSparseQRSupport.h.
{ eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()"); eigen_assert(b.cols()==1 && "This method is for vectors only"); //Compute Q^T * b typename Dest::PlainObject y, y2; y = matrixQ().transpose() * b; // Solves with the triangular matrix R Index rk = this->rank(); y2 = y; y.resize((std::max)(cols(),Index(y.rows())),y.cols()); y.topRows(rk) = this->matrixR().topLeftCorner(rk, rk).template triangularView<Upper>().solve(y2.topRows(rk)); // Apply the column permutation // colsPermutation() performs a copy of the permutation, // so let's apply it manually: for(Index i = 0; i < rk; ++i) dest.row(m_E[i]) = y.row(i); for(Index i = rk; i < cols(); ++i) dest.row(m_E[i]).setZero(); // y.bottomRows(y.rows()-rk).setZero(); // dest = colsPermutation() * y.topRows(cols()); m_info = Success; }
cholmod_common* Eigen::SPQR< _MatrixType >::cholmodCommon | ( | ) | const [inline] |
Definition at line 218 of file SuiteSparseQRSupport.h.
{ return &m_cc; }
Index Eigen::SPQR< _MatrixType >::cols | ( | void | ) | const [inline] |
Get the number of columns of the input matrix.
Definition at line 147 of file SuiteSparseQRSupport.h.
{ return m_cR->ncol; }
PermutationType Eigen::SPQR< _MatrixType >::colsPermutation | ( | ) | const [inline] |
Get the permutation that was applied to columns of A.
Definition at line 194 of file SuiteSparseQRSupport.h.
{ eigen_assert(m_isInitialized && "Decomposition is not initialized."); return PermutationType(m_E, m_cR->ncol); }
void Eigen::SPQR< _MatrixType >::compute | ( | const _MatrixType & | matrix | ) | [inline] |
Definition at line 103 of file SuiteSparseQRSupport.h.
{ if(m_isInitialized) SPQR_free(); MatrixType mat(matrix); /* Compute the default threshold as in MatLab, see: * Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing * Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011, Page 8:3 */ RealScalar pivotThreshold = m_tolerance; if(m_useDefaultThreshold) { RealScalar max2Norm = 0.0; for (int j = 0; j < mat.cols(); j++) max2Norm = numext::maxi(max2Norm, mat.col(j).norm()); if(max2Norm==RealScalar(0)) max2Norm = RealScalar(1); pivotThreshold = 20 * (mat.rows() + mat.cols()) * max2Norm * NumTraits<RealScalar>::epsilon(); } cholmod_sparse A; A = viewAsCholmod(mat); Index col = matrix.cols(); m_rank = SuiteSparseQR<Scalar>(m_ordering, pivotThreshold, col, &A, &m_cR, &m_E, &m_H, &m_HPinv, &m_HTau, &m_cc); if (!m_cR) { m_info = NumericalIssue; m_isInitialized = false; return; } m_info = Success; m_isInitialized = true; m_isRUpToDate = false; }
ComputationInfo Eigen::SPQR< _MatrixType >::info | ( | ) | const [inline] |
Reports whether previous computation was successful.
Success
if computation was succesful, NumericalIssue
if the sparse QR can not be computed Definition at line 226 of file SuiteSparseQRSupport.h.
{ eigen_assert(m_isInitialized && "Decomposition is not initialized."); return m_info; }
SPQRMatrixQReturnType<SPQR> Eigen::SPQR< _MatrixType >::matrixQ | ( | void | ) | const [inline] |
Get an expression of the matrix Q.
Definition at line 189 of file SuiteSparseQRSupport.h.
{
return SPQRMatrixQReturnType<SPQR>(*this);
}
const MatrixType Eigen::SPQR< _MatrixType >::matrixR | ( | ) | const [inline] |
Definition at line 179 of file SuiteSparseQRSupport.h.
{ eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()"); if(!m_isRUpToDate) { m_R = viewAsEigen<Scalar,ColMajor, typename MatrixType::StorageIndex>(*m_cR); m_isRUpToDate = true; } return m_R; }
Index Eigen::SPQR< _MatrixType >::rank | ( | ) | const [inline] |
Gets the rank of the matrix. It should be equal to matrixQR().cols if the matrix is full-rank
Definition at line 203 of file SuiteSparseQRSupport.h.
{ eigen_assert(m_isInitialized && "Decomposition is not initialized."); return m_cc.SPQR_istat[4]; }
Index Eigen::SPQR< _MatrixType >::rows | ( | void | ) | const [inline] |
Get the number of rows of the input matrix and the Q matrix
Definition at line 142 of file SuiteSparseQRSupport.h.
{return m_cR->nrow; }
void Eigen::SPQR< _MatrixType >::setPivotThreshold | ( | const RealScalar & | tol | ) | [inline] |
Set the tolerance tol to treat columns with 2-norm < =tol as zero.
Definition at line 211 of file SuiteSparseQRSupport.h.
{ m_useDefaultThreshold = false; m_tolerance = tol; }
void Eigen::SPQR< _MatrixType >::setSPQROrdering | ( | int | ord | ) | [inline] |
Set the fill-reducing ordering method to be used.
Definition at line 209 of file SuiteSparseQRSupport.h.
{ m_ordering = ord;}
void Eigen::SPQR< _MatrixType >::SPQR_free | ( | ) | [inline] |
friend struct SPQR_QProduct [friend] |
Definition at line 248 of file SuiteSparseQRSupport.h.
int Eigen::SPQR< _MatrixType >::m_allow_tol [protected] |
Definition at line 237 of file SuiteSparseQRSupport.h.
bool Eigen::SPQR< _MatrixType >::m_analysisIsOk [protected] |
Definition at line 232 of file SuiteSparseQRSupport.h.
cholmod_common Eigen::SPQR< _MatrixType >::m_cc [mutable, protected] |
Definition at line 246 of file SuiteSparseQRSupport.h.
cholmod_sparse* Eigen::SPQR< _MatrixType >::m_cR [mutable, protected] |
Definition at line 239 of file SuiteSparseQRSupport.h.
StorageIndex* Eigen::SPQR< _MatrixType >::m_E [mutable, protected] |
Definition at line 241 of file SuiteSparseQRSupport.h.
bool Eigen::SPQR< _MatrixType >::m_factorizationIsOk [protected] |
Definition at line 233 of file SuiteSparseQRSupport.h.
cholmod_sparse* Eigen::SPQR< _MatrixType >::m_H [mutable, protected] |
Definition at line 242 of file SuiteSparseQRSupport.h.
StorageIndex* Eigen::SPQR< _MatrixType >::m_HPinv [mutable, protected] |
Definition at line 243 of file SuiteSparseQRSupport.h.
cholmod_dense* Eigen::SPQR< _MatrixType >::m_HTau [mutable, protected] |
Definition at line 244 of file SuiteSparseQRSupport.h.
ComputationInfo Eigen::SPQR< _MatrixType >::m_info [mutable, protected] |
Definition at line 235 of file SuiteSparseQRSupport.h.
bool Eigen::SPQR< _MatrixType >::m_isRUpToDate [mutable, protected] |
Definition at line 234 of file SuiteSparseQRSupport.h.
int Eigen::SPQR< _MatrixType >::m_ordering [protected] |
Definition at line 236 of file SuiteSparseQRSupport.h.
MatrixType Eigen::SPQR< _MatrixType >::m_R [mutable, protected] |
Definition at line 240 of file SuiteSparseQRSupport.h.
Index Eigen::SPQR< _MatrixType >::m_rank [mutable, protected] |
Definition at line 245 of file SuiteSparseQRSupport.h.
RealScalar Eigen::SPQR< _MatrixType >::m_tolerance [protected] |
Definition at line 238 of file SuiteSparseQRSupport.h.
bool Eigen::SPQR< _MatrixType >::m_useDefaultThreshold [protected] |
Definition at line 247 of file SuiteSparseQRSupport.h.