Shapeworks Studio  2.1
Shape analysis software suite
List of all members | Public Types | Public Member Functions | Protected Attributes
Eigen::ColPivHouseholderQR< _MatrixType > Class Template Reference

Householder rank-revealing QR decomposition of a matrix with column-pivoting. More...

#include <ColPivHouseholderQR.h>

+ Collaboration diagram for Eigen::ColPivHouseholderQR< _MatrixType >:

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
}
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Index Index
 
typedef Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime > MatrixQType
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime > PermutationType
 
typedef internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
typedef internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
 
typedef HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::type > HouseholderSequenceType
 

Public Member Functions

 ColPivHouseholderQR ()
 Default Constructor. More...
 
 ColPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
 ColPivHouseholderQR (const MatrixType &matrix)
 Constructs a QR factorization from a given matrix. More...
 
template<typename Rhs >
const internal::solve_retval< ColPivHouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 
HouseholderSequenceType householderQ (void) const
 
HouseholderSequenceType matrixQ (void) const
 
const MatrixType & matrixQR () const
 
const MatrixType & matrixR () const
 
ColPivHouseholderQRcompute (const MatrixType &matrix)
 
const PermutationTypecolsPermutation () const
 
MatrixType::RealScalar absDeterminant () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
Index rank () const
 
Index dimensionOfKernel () const
 
bool isInjective () const
 
bool isSurjective () const
 
bool isInvertible () const
 
const internal::solve_retval< ColPivHouseholderQR, typename MatrixType::IdentityReturnType > inverse () const
 
Index rows () const
 
Index cols () const
 
const HCoeffsTypehCoeffs () const
 
ColPivHouseholderQRsetThreshold (const RealScalar &threshold)
 
ColPivHouseholderQRsetThreshold (Default_t)
 
RealScalar threshold () const
 
Index nonzeroPivots () const
 
RealScalar maxPivot () const
 
ComputationInfo info () const
 Reports whether the QR factorization was succesful. More...
 

Protected Attributes

MatrixType m_qr
 
HCoeffsType m_hCoeffs
 
PermutationType m_colsPermutation
 
IntRowVectorType m_colsTranspositions
 
RowVectorType m_temp
 
RealRowVectorType m_colSqNorms
 
bool m_isInitialized
 
bool m_usePrescribedThreshold
 
RealScalar m_prescribedThreshold
 
RealScalar m_maxpivot
 
Index m_nonzero_pivots
 
Index m_det_pq
 

Detailed Description

template<typename _MatrixType>
class Eigen::ColPivHouseholderQR< _MatrixType >

Householder rank-revealing QR decomposition of a matrix with column-pivoting.

Parameters
MatrixTypethe type of the matrix of which we are computing the QR decomposition

This class performs a rank-revealing QR decomposition of a matrix A into matrices P, Q and R such that

\[ \mathbf{A} \, \mathbf{P} = \mathbf{Q} \, \mathbf{R} \]

by using Householder transformations. Here, P is a permutation matrix, Q a unitary matrix and R an upper triangular matrix.

This decomposition performs column pivoting in order to be rank-revealing and improve numerical stability. It is slower than HouseholderQR, and faster than FullPivHouseholderQR.

See also
MatrixBase::colPivHouseholderQr()

Definition at line 222 of file ForwardDeclarations.h.

Constructor & Destructor Documentation

template<typename _MatrixType>
Eigen::ColPivHouseholderQR< _MatrixType >::ColPivHouseholderQR ( )
inline

Default Constructor.

The default constructor is useful in cases in which the user intends to perform decompositions via ColPivHouseholderQR::compute(const MatrixType&).

Definition at line 72 of file ColPivHouseholderQR.h.

73  : m_qr(),
74  m_hCoeffs(),
75  m_colsPermutation(),
76  m_colsTranspositions(),
77  m_temp(),
78  m_colSqNorms(),
79  m_isInitialized(false) {}
template<typename _MatrixType>
Eigen::ColPivHouseholderQR< _MatrixType >::ColPivHouseholderQR ( Index  rows,
Index  cols 
)
inline

Default Constructor with memory preallocation.

Like the default constructor but with preallocation of the internal data according to the specified problem size.

See also
ColPivHouseholderQR()

Definition at line 87 of file ColPivHouseholderQR.h.

88  : m_qr(rows, cols),
89  m_hCoeffs((std::min)(rows,cols)),
90  m_colsPermutation(PermIndexType(cols)),
91  m_colsTranspositions(cols),
92  m_temp(cols),
93  m_colSqNorms(cols),
94  m_isInitialized(false),
95  m_usePrescribedThreshold(false) {}
template<typename _MatrixType>
Eigen::ColPivHouseholderQR< _MatrixType >::ColPivHouseholderQR ( const MatrixType &  matrix)
inline

Constructs a QR factorization from a given matrix.

This constructor computes the QR factorization of the matrix matrix by calling the method compute(). It is a short cut for:

ColPivHouseholderQR<MatrixType> qr(matrix.rows(), matrix.cols());
qr.compute(matrix);
See also
compute()

Definition at line 109 of file ColPivHouseholderQR.h.

110  : m_qr(matrix.rows(), matrix.cols()),
111  m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
112  m_colsPermutation(PermIndexType(matrix.cols())),
113  m_colsTranspositions(matrix.cols()),
114  m_temp(matrix.cols()),
115  m_colSqNorms(matrix.cols()),
116  m_isInitialized(false),
117  m_usePrescribedThreshold(false)
118  {
119  compute(matrix);
120  }
Definition: math3d.h:219
ColPivHouseholderQR & compute(const MatrixType &matrix)

Member Function Documentation

template<typename MatrixType >
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType >::absDeterminant ( ) const
Returns
the absolute value of the determinant of the matrix of which *this is the QR decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the QR decomposition has already been computed.
Note
This is only for square matrices.
Warning
a determinant can be very big or small, so for matrices of large enough dimension, there is a risk of overflow/underflow. One way to work around that is to use logAbsDeterminant() instead.
See also
logAbsDeterminant(), MatrixBase::determinant()

Definition at line 399 of file ColPivHouseholderQR.h.

400 {
401  using std::abs;
402  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
403  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
404  return abs(m_qr.diagonal().prod());
405 }
template<typename _MatrixType>
const PermutationType& Eigen::ColPivHouseholderQR< _MatrixType >::colsPermutation ( ) const
inline
Returns
a const reference to the column permutation matrix

Definition at line 179 of file ColPivHouseholderQR.h.

180  {
181  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
182  return m_colsPermutation;
183  }
template<typename MatrixType >
ColPivHouseholderQR< MatrixType > & Eigen::ColPivHouseholderQR< MatrixType >::compute ( const MatrixType &  matrix)

Performs the QR factorization of the given matrix matrix. The result of the factorization is stored into *this, and a reference to *this is returned.

See also
class ColPivHouseholderQR, ColPivHouseholderQR(const MatrixType&)

Definition at line 422 of file ColPivHouseholderQR.h.

423 {
424  using std::abs;
425  Index rows = matrix.rows();
426  Index cols = matrix.cols();
427  Index size = matrix.diagonalSize();
428 
429  // the column permutation is stored as int indices, so just to be sure:
430  eigen_assert(cols<=NumTraits<int>::highest());
431 
432  m_qr = matrix;
433  m_hCoeffs.resize(size);
434 
435  m_temp.resize(cols);
436 
437  m_colsTranspositions.resize(matrix.cols());
438  Index number_of_transpositions = 0;
439 
440  m_colSqNorms.resize(cols);
441  for(Index k = 0; k < cols; ++k)
442  m_colSqNorms.coeffRef(k) = m_qr.col(k).squaredNorm();
443 
444  RealScalar threshold_helper = m_colSqNorms.maxCoeff() * numext::abs2(NumTraits<Scalar>::epsilon()) / RealScalar(rows);
445 
446  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
447  m_maxpivot = RealScalar(0);
448 
449  for(Index k = 0; k < size; ++k)
450  {
451  // first, we look up in our table m_colSqNorms which column has the biggest squared norm
452  Index biggest_col_index;
453  RealScalar biggest_col_sq_norm = m_colSqNorms.tail(cols-k).maxCoeff(&biggest_col_index);
454  biggest_col_index += k;
455 
456  // since our table m_colSqNorms accumulates imprecision at every step, we must now recompute
457  // the actual squared norm of the selected column.
458  // Note that not doing so does result in solve() sometimes returning inf/nan values
459  // when running the unit test with 1000 repetitions.
460  biggest_col_sq_norm = m_qr.col(biggest_col_index).tail(rows-k).squaredNorm();
461 
462  // we store that back into our table: it can't hurt to correct our table.
463  m_colSqNorms.coeffRef(biggest_col_index) = biggest_col_sq_norm;
464 
465  // if the current biggest column is smaller than epsilon times the initial biggest column,
466  // terminate to avoid generating nan/inf values.
467  // Note that here, if we test instead for "biggest == 0", we get a failure every 1000 (or so)
468  // repetitions of the unit test, with the result of solve() filled with large values of the order
469  // of 1/(size*epsilon).
470  if(biggest_col_sq_norm < threshold_helper * RealScalar(rows-k))
471  {
472  m_nonzero_pivots = k;
473  m_hCoeffs.tail(size-k).setZero();
474  m_qr.bottomRightCorner(rows-k,cols-k)
475  .template triangularView<StrictlyLower>()
476  .setZero();
477  break;
478  }
479 
480  // apply the transposition to the columns
481  m_colsTranspositions.coeffRef(k) = biggest_col_index;
482  if(k != biggest_col_index) {
483  m_qr.col(k).swap(m_qr.col(biggest_col_index));
484  std::swap(m_colSqNorms.coeffRef(k), m_colSqNorms.coeffRef(biggest_col_index));
485  ++number_of_transpositions;
486  }
487 
488  // generate the householder vector, store it below the diagonal
489  RealScalar beta;
490  m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
491 
492  // apply the householder transformation to the diagonal coefficient
493  m_qr.coeffRef(k,k) = beta;
494 
495  // remember the maximum absolute value of diagonal coefficients
496  if(abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
497 
498  // apply the householder transformation
499  m_qr.bottomRightCorner(rows-k, cols-k-1)
500  .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
501 
502  // update our table of squared norms of the columns
503  m_colSqNorms.tail(cols-k-1) -= m_qr.row(k).tail(cols-k-1).cwiseAbs2();
504  }
505 
506  m_colsPermutation.setIdentity(PermIndexType(cols));
507  for(PermIndexType k = 0; k < m_nonzero_pivots; ++k)
508  m_colsPermutation.applyTranspositionOnTheRight(k, PermIndexType(m_colsTranspositions.coeff(k)));
509 
510  m_det_pq = (number_of_transpositions%2) ? -1 : 1;
511  m_isInitialized = true;
512 
513  return *this;
514 }
Definition: math3d.h:219
Derived & applyTranspositionOnTheRight(Index i, Index j)
template<typename _MatrixType>
Index Eigen::ColPivHouseholderQR< _MatrixType >::dimensionOfKernel ( ) const
inline
Returns
the dimension of the kernel of the matrix of which *this is the QR decomposition.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 237 of file ColPivHouseholderQR.h.

238  {
239  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
240  return cols() - rank();
241  }
template<typename _MatrixType>
const HCoeffsType& Eigen::ColPivHouseholderQR< _MatrixType >::hCoeffs ( ) const
inline
Returns
a const reference to the vector of Householder coefficients used to represent the factor Q.

For advanced uses only.

Definition at line 302 of file ColPivHouseholderQR.h.

302 { return m_hCoeffs; }
template<typename MatrixType >
ColPivHouseholderQR< MatrixType >::HouseholderSequenceType Eigen::ColPivHouseholderQR< MatrixType >::householderQ ( void  ) const
Returns
the matrix Q as a sequence of householder transformations

Definition at line 560 of file ColPivHouseholderQR.h.

561 {
562  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
563  return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate()).setLength(m_nonzero_pivots);
564 }
template<typename _MatrixType>
ComputationInfo Eigen::ColPivHouseholderQR< _MatrixType >::info ( ) const
inline

Reports whether the QR factorization was succesful.

Note
This function always returns Success. It is provided for compatibility with other factorization routines.
Returns
Success

Definition at line 379 of file ColPivHouseholderQR.h.

380  {
381  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
382  return Success;
383  }
template<typename _MatrixType>
const internal::solve_retval<ColPivHouseholderQR, typename MatrixType::IdentityReturnType> Eigen::ColPivHouseholderQR< _MatrixType >::inverse ( void  ) const
inline
Returns
the inverse of the matrix of which *this is the QR decomposition.
Note
If this matrix is not invertible, the returned matrix has undefined coefficients. Use isInvertible() to first determine whether this matrix is invertible.

Definition at line 288 of file ColPivHouseholderQR.h.

289  {
290  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
291  return internal::solve_retval<ColPivHouseholderQR,typename MatrixType::IdentityReturnType>
292  (*this, MatrixType::Identity(m_qr.rows(), m_qr.cols()));
293  }
template<typename _MatrixType>
bool Eigen::ColPivHouseholderQR< _MatrixType >::isInjective ( ) const
inline
Returns
true if the matrix of which *this is the QR decomposition represents an injective linear map, i.e. has trivial kernel; false otherwise.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 250 of file ColPivHouseholderQR.h.

251  {
252  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
253  return rank() == cols();
254  }
template<typename _MatrixType>
bool Eigen::ColPivHouseholderQR< _MatrixType >::isInvertible ( ) const
inline
Returns
true if the matrix of which *this is the QR decomposition is invertible.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 275 of file ColPivHouseholderQR.h.

276  {
277  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
278  return isInjective() && isSurjective();
279  }
template<typename _MatrixType>
bool Eigen::ColPivHouseholderQR< _MatrixType >::isSurjective ( ) const
inline
Returns
true if the matrix of which *this is the QR decomposition represents a surjective linear map; false otherwise.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 263 of file ColPivHouseholderQR.h.

264  {
265  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
266  return rank() == rows();
267  }
template<typename MatrixType >
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType >::logAbsDeterminant ( ) const
Returns
the natural log of the absolute value of the determinant of the matrix of which *this is the QR decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the QR decomposition has already been computed.
Note
This is only for square matrices.
This method is useful to work around the risk of overflow/underflow that's inherent to determinant computation.
See also
absDeterminant(), MatrixBase::determinant()

Definition at line 408 of file ColPivHouseholderQR.h.

409 {
410  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
411  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
412  return m_qr.diagonal().cwiseAbs().array().log().sum();
413 }
template<typename _MatrixType>
const MatrixType& Eigen::ColPivHouseholderQR< _MatrixType >::matrixQR ( ) const
inline
Returns
a reference to the matrix where the Householder QR decomposition is stored

Definition at line 155 of file ColPivHouseholderQR.h.

156  {
157  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
158  return m_qr;
159  }
template<typename _MatrixType>
const MatrixType& Eigen::ColPivHouseholderQR< _MatrixType >::matrixR ( void  ) const
inline
Returns
a reference to the matrix where the result Householder QR is stored
Warning
The strict lower part of this matrix contains internal values. Only the upper triangular part should be referenced. To get it, use
matrixR().template triangularView<Upper>()
For rank-deficient matrices, use
matrixR().topLeftCorner(rank(), rank()).template triangularView<Upper>()

Definition at line 170 of file ColPivHouseholderQR.h.

171  {
172  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
173  return m_qr;
174  }
template<typename _MatrixType>
RealScalar Eigen::ColPivHouseholderQR< _MatrixType >::maxPivot ( ) const
inline
Returns
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of R.

Definition at line 371 of file ColPivHouseholderQR.h.

371 { return m_maxpivot; }
template<typename _MatrixType>
Index Eigen::ColPivHouseholderQR< _MatrixType >::nonzeroPivots ( ) const
inline
Returns
the number of nonzero pivots in the QR decomposition. Here nonzero is meant in the exact sense, not in a fuzzy sense. So that notion isn't really intrinsically interesting, but it is still useful when implementing algorithms.
See also
rank()

Definition at line 362 of file ColPivHouseholderQR.h.

363  {
364  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
365  return m_nonzero_pivots;
366  }
template<typename _MatrixType>
Index Eigen::ColPivHouseholderQR< _MatrixType >::rank ( ) const
inline
Returns
the rank of the matrix of which *this is the QR decomposition.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 220 of file ColPivHouseholderQR.h.

221  {
222  using std::abs;
223  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
224  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
225  Index result = 0;
226  for(Index i = 0; i < m_nonzero_pivots; ++i)
227  result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold);
228  return result;
229  }
template<typename _MatrixType>
ColPivHouseholderQR& Eigen::ColPivHouseholderQR< _MatrixType >::setThreshold ( const RealScalar &  threshold)
inline

Allows to prescribe a threshold to be used by certain methods, such as rank(), who need to determine when pivots are to be considered nonzero. This is not used for the QR decomposition itself.

When it needs to get the threshold value, Eigen calls threshold(). By default, this uses a formula to automatically determine a reasonable threshold. Once you have called the present method setThreshold(const RealScalar&), your value is used instead.

Parameters
thresholdThe new value to use as the threshold.

A pivot will be considered nonzero if its absolute value is strictly greater than $ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert $ where maxpivot is the biggest pivot.

If you want to come back to the default behavior, call setThreshold(Default_t)

Definition at line 321 of file ColPivHouseholderQR.h.

322  {
323  m_usePrescribedThreshold = true;
324  m_prescribedThreshold = threshold;
325  return *this;
326  }
template<typename _MatrixType>
ColPivHouseholderQR& Eigen::ColPivHouseholderQR< _MatrixType >::setThreshold ( Default_t  )
inline

Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold.

You should pass the special object Eigen::Default as parameter here.

qr.setThreshold(Eigen::Default);

See the documentation of setThreshold(const RealScalar&).

Definition at line 336 of file ColPivHouseholderQR.h.

337  {
338  m_usePrescribedThreshold = false;
339  return *this;
340  }
template<typename _MatrixType>
template<typename Rhs >
const internal::solve_retval<ColPivHouseholderQR, Rhs> Eigen::ColPivHouseholderQR< _MatrixType >::solve ( const MatrixBase< Rhs > &  b) const
inline

This method finds a solution x to the equation Ax=b, where A is the matrix of which *this is the QR decomposition, if any exists.

Parameters
bthe right-hand-side of the equation to solve.
Returns
a solution.
Note
The case where b is a matrix is not yet implemented. Also, this code is space inefficient.

Example:

Output:

 

Definition at line 141 of file ColPivHouseholderQR.h.

142  {
143  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
144  return internal::solve_retval<ColPivHouseholderQR, Rhs>(*this, b.derived());
145  }
template<typename _MatrixType>
RealScalar Eigen::ColPivHouseholderQR< _MatrixType >::threshold ( ) const
inline

Returns the threshold that will be used by certain methods such as rank().

See the documentation of setThreshold(const RealScalar&).

Definition at line 346 of file ColPivHouseholderQR.h.

347  {
348  eigen_assert(m_isInitialized || m_usePrescribedThreshold);
349  return m_usePrescribedThreshold ? m_prescribedThreshold
350  // this formula comes from experimenting (see "LU precision tuning" thread on the list)
351  // and turns out to be identical to Higham's formula used already in LDLt.
352  : NumTraits<Scalar>::epsilon() * m_qr.diagonalSize();
353  }

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