Eigen::SPQR< MatrixType_ > Class Template Reference

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

+ Inheritance diagram for Eigen::SPQR< MatrixType_ >:

Public Types

enum  {
  ColsAtCompileTime ,
  MaxColsAtCompileTime
}
 
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 ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< SPQR< MatrixType_ > >
SPQR< MatrixType_ > & derived ()
 
const SPQR< MatrixType_ > & derived () const
 
const Solve< SPQR< MatrixType_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< SPQR< MatrixType_ >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

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_isRUpToDate
 
int m_ordering
 
MatrixType m_R
 
Index m_rank
 
Index m_rows
 
RealScalar m_tolerance
 
bool m_useDefaultThreshold
 
- Protected Attributes inherited from Eigen::SparseSolverBase< SPQR< MatrixType_ > >
bool m_isInitialized
 

Detailed Description

template<typename MatrixType_>
class Eigen::SPQR< MatrixType_ >

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
MatrixType_The type of the sparse matrix A, must be a column-major SparseMatrix<>

This class follows the sparse solver concept .

Definition at line 62 of file SuiteSparseQRSupport.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ >
typedef SparseSolverBase<SPQR<MatrixType_> > Eigen::SPQR< MatrixType_ >::Base
protected

Definition at line 65 of file SuiteSparseQRSupport.h.

◆ MatrixType

template<typename MatrixType_ >
typedef SparseMatrix<Scalar, ColMajor, StorageIndex> Eigen::SPQR< MatrixType_ >::MatrixType

Definition at line 71 of file SuiteSparseQRSupport.h.

◆ PermutationType

template<typename MatrixType_ >
typedef Map<PermutationMatrix<Dynamic, Dynamic, StorageIndex> > Eigen::SPQR< MatrixType_ >::PermutationType

Definition at line 72 of file SuiteSparseQRSupport.h.

◆ RealScalar

template<typename MatrixType_ >
typedef MatrixType_::RealScalar Eigen::SPQR< MatrixType_ >::RealScalar

Definition at line 69 of file SuiteSparseQRSupport.h.

◆ Scalar

template<typename MatrixType_ >
typedef MatrixType_::Scalar Eigen::SPQR< MatrixType_ >::Scalar

Definition at line 68 of file SuiteSparseQRSupport.h.

◆ StorageIndex

template<typename MatrixType_ >
typedef SuiteSparse_long Eigen::SPQR< MatrixType_ >::StorageIndex

Definition at line 70 of file SuiteSparseQRSupport.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 73 of file SuiteSparseQRSupport.h.

73  {
76  };
const int Dynamic
Definition: Constants.h:24

Constructor & Destructor Documentation

◆ SPQR() [1/2]

template<typename MatrixType_ >
Eigen::SPQR< MatrixType_ >::SPQR ( )
inline

Definition at line 78 of file SuiteSparseQRSupport.h.

79  : m_analysisIsOk(false),
80  m_factorizationIsOk(false),
81  m_isRUpToDate(false),
82  m_ordering(SPQR_ORDERING_DEFAULT),
83  m_allow_tol(SPQR_DEFAULT_TOL),
84  m_tolerance (NumTraits<Scalar>::epsilon()),
85  m_cR(0),
86  m_E(0),
87  m_H(0),
88  m_HPinv(0),
89  m_HTau(0),
91  {
92  cholmod_l_start(&m_cc);
93  }
cholmod_dense * m_HTau
StorageIndex * m_E
cholmod_sparse * m_cR
RealScalar m_tolerance
cholmod_sparse * m_H
StorageIndex * m_HPinv
cholmod_common m_cc

◆ SPQR() [2/2]

template<typename MatrixType_ >
Eigen::SPQR< MatrixType_ >::SPQR ( const MatrixType_ &  matrix)
inlineexplicit

Definition at line 95 of file SuiteSparseQRSupport.h.

96  : m_analysisIsOk(false),
97  m_factorizationIsOk(false),
98  m_isRUpToDate(false),
99  m_ordering(SPQR_ORDERING_DEFAULT),
100  m_allow_tol(SPQR_DEFAULT_TOL),
101  m_tolerance (NumTraits<Scalar>::epsilon()),
102  m_cR(0),
103  m_E(0),
104  m_H(0),
105  m_HPinv(0),
106  m_HTau(0),
108  {
109  cholmod_l_start(&m_cc);
110  compute(matrix);
111  }
void compute(const MatrixType_ &matrix)

◆ ~SPQR()

template<typename MatrixType_ >
Eigen::SPQR< MatrixType_ >::~SPQR ( )
inline

Definition at line 113 of file SuiteSparseQRSupport.h.

114  {
115  SPQR_free();
116  cholmod_l_finish(&m_cc);
117  }

Member Function Documentation

◆ _solve_impl()

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

Definition at line 174 of file SuiteSparseQRSupport.h.

175  {
176  eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
177  eigen_assert(b.cols()==1 && "This method is for vectors only");
178 
179  //Compute Q^T * b
180  typename Dest::PlainObject y, y2;
181  y = matrixQ().transpose() * b;
182 
183  // Solves with the triangular matrix R
184  Index rk = this->rank();
185  y2 = y;
186  y.resize((std::max)(cols(),Index(y.rows())),y.cols());
187  y.topRows(rk) = this->matrixR().topLeftCorner(rk, rk).template triangularView<Upper>().solve(y2.topRows(rk));
188 
189  // Apply the column permutation
190  // colsPermutation() performs a copy of the permutation,
191  // so let's apply it manually:
192  for(Index i = 0; i < rk; ++i) dest.row(m_E[i]) = y.row(i);
193  for(Index i = rk; i < cols(); ++i) dest.row(m_E[i]).setZero();
194 
195 // y.bottomRows(y.rows()-rk).setZero();
196 // dest = colsPermutation() * y.topRows(cols());
197 
198  m_info = Success;
199  }
Array< int, 3, 1 > b
#define eigen_assert(x)
Definition: Macros.h:902
Index rank() const
const MatrixType matrixR() const
ComputationInfo m_info
Index cols() const
SPQRMatrixQReturnType< SPQR > matrixQ() const
Get an expression of the matrix Q.
@ Success
Definition: Constants.h:446
bfloat16() max(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:690
const Scalar & y
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82

◆ cholmodCommon()

template<typename MatrixType_ >
cholmod_common* Eigen::SPQR< MatrixType_ >::cholmodCommon ( ) const
inline
Returns
a pointer to the SPQR workspace

Definition at line 242 of file SuiteSparseQRSupport.h.

242 { return &m_cc; }

◆ cols()

template<typename MatrixType_ >
Index Eigen::SPQR< MatrixType_ >::cols ( void  ) const
inline

Get the number of columns of the input matrix.

Definition at line 171 of file SuiteSparseQRSupport.h.

171 { return m_cR->ncol; }

◆ colsPermutation()

template<typename MatrixType_ >
PermutationType Eigen::SPQR< MatrixType_ >::colsPermutation ( ) const
inline

Get the permutation that was applied to columns of A.

Definition at line 218 of file SuiteSparseQRSupport.h.

219  {
220  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
221  return PermutationType(m_E, m_cR->ncol);
222  }
Map< PermutationMatrix< Dynamic, Dynamic, StorageIndex > > PermutationType

◆ compute()

template<typename MatrixType_ >
void Eigen::SPQR< MatrixType_ >::compute ( const MatrixType_ &  matrix)
inline

Definition at line 127 of file SuiteSparseQRSupport.h.

128  {
130 
131  MatrixType mat(matrix);
132 
133  /* Compute the default threshold as in MatLab, see:
134  * Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing
135  * Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011, Page 8:3
136  */
137  RealScalar pivotThreshold = m_tolerance;
139  {
140  RealScalar max2Norm = 0.0;
141  for (int j = 0; j < mat.cols(); j++) max2Norm = numext::maxi(max2Norm, mat.col(j).norm());
142  if(numext::is_exactly_zero(max2Norm))
143  max2Norm = RealScalar(1);
144  pivotThreshold = 20 * (mat.rows() + mat.cols()) * max2Norm * NumTraits<RealScalar>::epsilon();
145  }
146  cholmod_sparse A;
147  A = viewAsCholmod(mat);
148  m_rows = matrix.rows();
149  Index col = matrix.cols();
150  m_rank = SuiteSparseQR<Scalar>(m_ordering, pivotThreshold, col, &A,
151  &m_cR, &m_E, &m_H, &m_HPinv, &m_HTau, &m_cc);
152 
153  if (!m_cR)
154  {
156  m_isInitialized = false;
157  return;
158  }
159  m_info = Success;
160  m_isInitialized = true;
161  m_isRUpToDate = false;
162  }
ColXpr col(Index i)
This is the const version of col().
MatrixXcf A
Matrix< float, 1, Dynamic > MatrixType
MatrixType_::RealScalar RealScalar
@ NumericalIssue
Definition: Constants.h:448
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< Scalar_, Options_, StorageIndex_ > > mat)
std::ptrdiff_t j

◆ info()

template<typename MatrixType_ >
ComputationInfo Eigen::SPQR< MatrixType_ >::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 250 of file SuiteSparseQRSupport.h.

251  {
252  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
253  return m_info;
254  }

◆ matrixQ()

template<typename MatrixType_ >
SPQRMatrixQReturnType<SPQR> Eigen::SPQR< MatrixType_ >::matrixQ ( void  ) const
inline

Get an expression of the matrix Q.

Definition at line 213 of file SuiteSparseQRSupport.h.

214  {
215  return SPQRMatrixQReturnType<SPQR>(*this);
216  }

◆ matrixR()

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

Definition at line 203 of file SuiteSparseQRSupport.h.

204  {
205  eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
206  if(!m_isRUpToDate) {
207  m_R = viewAsEigen<Scalar,ColMajor, typename MatrixType::StorageIndex>(*m_cR);
208  m_isRUpToDate = true;
209  }
210  return m_R;
211  }

◆ rank()

template<typename MatrixType_ >
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 227 of file SuiteSparseQRSupport.h.

228  {
229  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
230  return m_cc.SPQR_istat[4];
231  }

◆ rows()

template<typename MatrixType_ >
Index Eigen::SPQR< MatrixType_ >::rows ( void  ) const
inline

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

Definition at line 166 of file SuiteSparseQRSupport.h.

166 {return m_rows; }

◆ setPivotThreshold()

template<typename MatrixType_ >
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 235 of file SuiteSparseQRSupport.h.

236  {
237  m_useDefaultThreshold = false;
238  m_tolerance = tol;
239  }

◆ setSPQROrdering()

template<typename MatrixType_ >
void Eigen::SPQR< MatrixType_ >::setSPQROrdering ( int  ord)
inline

Set the fill-reducing ordering method to be used.

Definition at line 233 of file SuiteSparseQRSupport.h.

233 { m_ordering = ord;}

◆ SPQR_free()

template<typename MatrixType_ >
void Eigen::SPQR< MatrixType_ >::SPQR_free ( )
inline

Definition at line 118 of file SuiteSparseQRSupport.h.

119  {
120  cholmod_l_free_sparse(&m_H, &m_cc);
121  cholmod_l_free_sparse(&m_cR, &m_cc);
122  cholmod_l_free_dense(&m_HTau, &m_cc);
123  std::free(m_E);
124  std::free(m_HPinv);
125  }

Member Data Documentation

◆ m_allow_tol

template<typename MatrixType_ >
int Eigen::SPQR< MatrixType_ >::m_allow_tol
protected

Definition at line 261 of file SuiteSparseQRSupport.h.

◆ m_analysisIsOk

template<typename MatrixType_ >
bool Eigen::SPQR< MatrixType_ >::m_analysisIsOk
protected

Definition at line 256 of file SuiteSparseQRSupport.h.

◆ m_cc

template<typename MatrixType_ >
cholmod_common Eigen::SPQR< MatrixType_ >::m_cc
mutableprotected

Definition at line 270 of file SuiteSparseQRSupport.h.

◆ m_cR

template<typename MatrixType_ >
cholmod_sparse* Eigen::SPQR< MatrixType_ >::m_cR
mutableprotected

Definition at line 263 of file SuiteSparseQRSupport.h.

◆ m_E

template<typename MatrixType_ >
StorageIndex* Eigen::SPQR< MatrixType_ >::m_E
mutableprotected

Definition at line 265 of file SuiteSparseQRSupport.h.

◆ m_factorizationIsOk

template<typename MatrixType_ >
bool Eigen::SPQR< MatrixType_ >::m_factorizationIsOk
protected

Definition at line 257 of file SuiteSparseQRSupport.h.

◆ m_H

template<typename MatrixType_ >
cholmod_sparse* Eigen::SPQR< MatrixType_ >::m_H
mutableprotected

Definition at line 266 of file SuiteSparseQRSupport.h.

◆ m_HPinv

template<typename MatrixType_ >
StorageIndex* Eigen::SPQR< MatrixType_ >::m_HPinv
mutableprotected

Definition at line 267 of file SuiteSparseQRSupport.h.

◆ m_HTau

template<typename MatrixType_ >
cholmod_dense* Eigen::SPQR< MatrixType_ >::m_HTau
mutableprotected

Definition at line 268 of file SuiteSparseQRSupport.h.

◆ m_info

template<typename MatrixType_ >
ComputationInfo Eigen::SPQR< MatrixType_ >::m_info
mutableprotected

Definition at line 259 of file SuiteSparseQRSupport.h.

◆ m_isRUpToDate

template<typename MatrixType_ >
bool Eigen::SPQR< MatrixType_ >::m_isRUpToDate
mutableprotected

Definition at line 258 of file SuiteSparseQRSupport.h.

◆ m_ordering

template<typename MatrixType_ >
int Eigen::SPQR< MatrixType_ >::m_ordering
protected

Definition at line 260 of file SuiteSparseQRSupport.h.

◆ m_R

template<typename MatrixType_ >
MatrixType Eigen::SPQR< MatrixType_ >::m_R
mutableprotected

Definition at line 264 of file SuiteSparseQRSupport.h.

◆ m_rank

template<typename MatrixType_ >
Index Eigen::SPQR< MatrixType_ >::m_rank
mutableprotected

Definition at line 269 of file SuiteSparseQRSupport.h.

◆ m_rows

template<typename MatrixType_ >
Index Eigen::SPQR< MatrixType_ >::m_rows
protected

Definition at line 272 of file SuiteSparseQRSupport.h.

◆ m_tolerance

template<typename MatrixType_ >
RealScalar Eigen::SPQR< MatrixType_ >::m_tolerance
protected

Definition at line 262 of file SuiteSparseQRSupport.h.

◆ m_useDefaultThreshold

template<typename MatrixType_ >
bool Eigen::SPQR< MatrixType_ >::m_useDefaultThreshold
protected

Definition at line 271 of file SuiteSparseQRSupport.h.


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