Public Types | Public Member Functions | Protected Types | Protected Attributes | Friends | List of all members
Eigen::SPQR Class Reference

Sparse QR factorization based on SuiteSparseQR library. More...

#include <SuiteSparseQRSupport.h>

Public Types

enum  { ColsAtCompileTime = Dynamic, MaxColsAtCompileTime = Dynamic }
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexMatrixType
 
typedef Map< PermutationMatrix< Dynamic, Dynamic, StorageIndex > > PermutationType
 
typedef _MatrixType::RealScalar RealScalar
 
typedef _MatrixType::Scalar Scalar
 
typedef SuiteSparse_long StorageIndex
 

Public Member Functions

template<typename Rhs , typename Dest >
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
cholmod_common * cholmodCommon () const
 
Index cols () const
 
PermutationType colsPermutation () const
 Get the permutation that was applied to columns of A. More...
 
void compute (const _MatrixType &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
SPQRMatrixQReturnType< SPQRmatrixQ () const
 Get an expression of the matrix Q. More...
 
const MatrixType matrixR () const
 
Index rank () const
 
Index rows () const
 
void setPivotThreshold (const RealScalar &tol)
 Set the tolerance tol to treat columns with 2-norm < =tol as zero. More...
 
void setSPQROrdering (int ord)
 Set the fill-reducing ordering method to be used. More...
 
 SPQR ()
 
 SPQR (const _MatrixType &matrix)
 
void SPQR_free ()
 
 ~SPQR ()
 

Protected Types

typedef SparseSolverBase< SPQR< _MatrixType > > Base
 

Protected Attributes

int m_allow_tol
 
bool m_analysisIsOk
 
cholmod_common m_cc
 
cholmod_sparse * m_cR
 
StorageIndexm_E
 
bool m_factorizationIsOk
 
cholmod_sparse * m_H
 
StorageIndexm_HPinv
 
cholmod_dense * m_HTau
 
ComputationInfo m_info
 
bool m_isInitialized
 
bool m_isRUpToDate
 
int m_ordering
 
MatrixType m_R
 
Index m_rank
 
Index m_rows
 
RealScalar m_tolerance
 
bool m_useDefaultThreshold
 

Friends

template<typename , typename >
struct SPQR_QProduct
 

Detailed Description

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

Template Parameters
_MatrixTypeThe type of the sparse matrix A, must be a column-major SparseMatrix<>

\implsparsesolverconcept

Definition at line 16 of file SuiteSparseQRSupport.h.

Member Typedef Documentation

◆ Base

typedef SparseSolverBase<SPQR<_MatrixType> > Eigen::SPQR::Base
protected

Definition at line 63 of file SuiteSparseQRSupport.h.

◆ MatrixType

Definition at line 69 of file SuiteSparseQRSupport.h.

◆ PermutationType

Definition at line 70 of file SuiteSparseQRSupport.h.

◆ RealScalar

typedef _MatrixType::RealScalar Eigen::SPQR::RealScalar

Definition at line 67 of file SuiteSparseQRSupport.h.

◆ Scalar

typedef _MatrixType::Scalar Eigen::SPQR::Scalar

Definition at line 66 of file SuiteSparseQRSupport.h.

◆ StorageIndex

Definition at line 68 of file SuiteSparseQRSupport.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 71 of file SuiteSparseQRSupport.h.

Constructor & Destructor Documentation

◆ SPQR() [1/2]

Eigen::SPQR::SPQR ( )
inline

Definition at line 76 of file SuiteSparseQRSupport.h.

◆ SPQR() [2/2]

Eigen::SPQR::SPQR ( const _MatrixType &  matrix)
inlineexplicit

Definition at line 93 of file SuiteSparseQRSupport.h.

◆ ~SPQR()

Eigen::SPQR::~SPQR ( )
inline

Definition at line 111 of file SuiteSparseQRSupport.h.

Member Function Documentation

◆ _solve_impl()

template<typename Rhs , typename Dest >
void Eigen::SPQR::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  dest 
) const
inline

Definition at line 172 of file SuiteSparseQRSupport.h.

◆ cholmodCommon()

cholmod_common* Eigen::SPQR::cholmodCommon ( ) const
inline
Returns
a pointer to the SPQR workspace

Definition at line 240 of file SuiteSparseQRSupport.h.

◆ cols()

Index Eigen::SPQR::cols ( ) const
inline

Get the number of columns of the input matrix.

Definition at line 169 of file SuiteSparseQRSupport.h.

◆ colsPermutation()

PermutationType Eigen::SPQR::colsPermutation ( ) const
inline

Get the permutation that was applied to columns of A.

Definition at line 216 of file SuiteSparseQRSupport.h.

◆ compute()

void Eigen::SPQR::compute ( const _MatrixType &  matrix)
inline

Definition at line 125 of file SuiteSparseQRSupport.h.

◆ info()

ComputationInfo Eigen::SPQR::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the sparse QR can not be computed

Definition at line 248 of file SuiteSparseQRSupport.h.

◆ matrixQ()

SPQRMatrixQReturnType<SPQR> Eigen::SPQR::matrixQ ( ) const
inline

Get an expression of the matrix Q.

Definition at line 211 of file SuiteSparseQRSupport.h.

◆ matrixR()

const MatrixType Eigen::SPQR::matrixR ( ) const
inline
Returns
the sparse triangular factor R. It is a sparse matrix

Definition at line 201 of file SuiteSparseQRSupport.h.

◆ rank()

Index Eigen::SPQR::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 225 of file SuiteSparseQRSupport.h.

◆ rows()

Index Eigen::SPQR::rows ( ) const
inline

Get the number of rows of the input matrix and the Q matrix

Definition at line 164 of file SuiteSparseQRSupport.h.

◆ setPivotThreshold()

void Eigen::SPQR::setPivotThreshold ( const RealScalar tol)
inline

Set the tolerance tol to treat columns with 2-norm < =tol as zero.

Definition at line 233 of file SuiteSparseQRSupport.h.

◆ setSPQROrdering()

void Eigen::SPQR::setSPQROrdering ( int  ord)
inline

Set the fill-reducing ordering method to be used.

Definition at line 231 of file SuiteSparseQRSupport.h.

◆ SPQR_free()

void Eigen::SPQR::SPQR_free ( )
inline

Definition at line 116 of file SuiteSparseQRSupport.h.

Friends And Related Function Documentation

◆ SPQR_QProduct

template<typename , typename >
friend struct SPQR_QProduct
friend

Definition at line 271 of file SuiteSparseQRSupport.h.

Member Data Documentation

◆ m_allow_tol

int Eigen::SPQR::m_allow_tol
protected

Definition at line 259 of file SuiteSparseQRSupport.h.

◆ m_analysisIsOk

bool Eigen::SPQR::m_analysisIsOk
protected

Definition at line 254 of file SuiteSparseQRSupport.h.

◆ m_cc

cholmod_common Eigen::SPQR::m_cc
mutableprotected

Definition at line 268 of file SuiteSparseQRSupport.h.

◆ m_cR

cholmod_sparse* Eigen::SPQR::m_cR
mutableprotected

Definition at line 261 of file SuiteSparseQRSupport.h.

◆ m_E

StorageIndex* Eigen::SPQR::m_E
mutableprotected

Definition at line 263 of file SuiteSparseQRSupport.h.

◆ m_factorizationIsOk

bool Eigen::SPQR::m_factorizationIsOk
protected

Definition at line 255 of file SuiteSparseQRSupport.h.

◆ m_H

cholmod_sparse* Eigen::SPQR::m_H
mutableprotected

Definition at line 264 of file SuiteSparseQRSupport.h.

◆ m_HPinv

StorageIndex* Eigen::SPQR::m_HPinv
mutableprotected

Definition at line 265 of file SuiteSparseQRSupport.h.

◆ m_HTau

cholmod_dense* Eigen::SPQR::m_HTau
mutableprotected

Definition at line 266 of file SuiteSparseQRSupport.h.

◆ m_info

ComputationInfo Eigen::SPQR::m_info
mutableprotected

Definition at line 257 of file SuiteSparseQRSupport.h.

◆ m_isInitialized

bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

Definition at line 119 of file SparseSolverBase.h.

◆ m_isRUpToDate

bool Eigen::SPQR::m_isRUpToDate
mutableprotected

Definition at line 256 of file SuiteSparseQRSupport.h.

◆ m_ordering

int Eigen::SPQR::m_ordering
protected

Definition at line 258 of file SuiteSparseQRSupport.h.

◆ m_R

MatrixType Eigen::SPQR::m_R
mutableprotected

Definition at line 262 of file SuiteSparseQRSupport.h.

◆ m_rank

Index Eigen::SPQR::m_rank
mutableprotected

Definition at line 267 of file SuiteSparseQRSupport.h.

◆ m_rows

Index Eigen::SPQR::m_rows
protected

Definition at line 270 of file SuiteSparseQRSupport.h.

◆ m_tolerance

RealScalar Eigen::SPQR::m_tolerance
protected

Definition at line 260 of file SuiteSparseQRSupport.h.

◆ m_useDefaultThreshold

bool Eigen::SPQR::m_useDefaultThreshold
protected

Definition at line 269 of file SuiteSparseQRSupport.h.


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


gtsam
Author(s):
autogenerated on Wed May 15 2024 15:29:28