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

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

#include <FullPivHouseholderQR.h>

+ Collaboration diagram for Eigen::FullPivHouseholderQR< _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 internal::FullPivHouseholderQRMatrixQReturnType< MatrixType > MatrixQReturnType
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef Matrix< Index, 1, ColsAtCompileTime, RowMajor, 1, MaxColsAtCompileTime > IntRowVectorType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime > PermutationType
 
typedef internal::plain_col_type< MatrixType, Index >::type IntColVectorType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
typedef internal::plain_col_type< MatrixType >::type ColVectorType
 

Public Member Functions

 FullPivHouseholderQR ()
 Default Constructor. More...
 
 FullPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
 FullPivHouseholderQR (const MatrixType &matrix)
 Constructs a QR factorization from a given matrix. More...
 
template<typename Rhs >
const internal::solve_retval< FullPivHouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 
MatrixQReturnType matrixQ (void) const
 
const MatrixType & matrixQR () const
 
FullPivHouseholderQRcompute (const MatrixType &matrix)
 
const PermutationTypecolsPermutation () const
 
const IntColVectorTyperowsTranspositions () 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< FullPivHouseholderQR, typename MatrixType::IdentityReturnType > inverse () const
 
Index rows () const
 
Index cols () const
 
const HCoeffsTypehCoeffs () const
 
FullPivHouseholderQRsetThreshold (const RealScalar &threshold)
 
FullPivHouseholderQRsetThreshold (Default_t)
 
RealScalar threshold () const
 
Index nonzeroPivots () const
 
RealScalar maxPivot () const
 

Protected Attributes

MatrixType m_qr
 
HCoeffsType m_hCoeffs
 
IntColVectorType m_rows_transpositions
 
IntRowVectorType m_cols_transpositions
 
PermutationType m_cols_permutation
 
RowVectorType m_temp
 
bool m_isInitialized
 
bool m_usePrescribedThreshold
 
RealScalar m_prescribedThreshold
 
RealScalar m_maxpivot
 
Index m_nonzero_pivots
 
RealScalar m_precision
 
Index m_det_pq
 

Detailed Description

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

Householder rank-revealing QR decomposition of a matrix with full 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 a very prudent full pivoting in order to be rank-revealing and achieve optimal numerical stability. The trade-off is that it is slower than HouseholderQR and ColPivHouseholderQR.

See also
MatrixBase::fullPivHouseholderQr()

Definition at line 223 of file ForwardDeclarations.h.

Constructor & Destructor Documentation

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

Default Constructor.

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

Definition at line 77 of file FullPivHouseholderQR.h.

78  : m_qr(),
79  m_hCoeffs(),
80  m_rows_transpositions(),
81  m_cols_transpositions(),
82  m_cols_permutation(),
83  m_temp(),
84  m_isInitialized(false),
85  m_usePrescribedThreshold(false) {}
template<typename _MatrixType>
Eigen::FullPivHouseholderQR< _MatrixType >::FullPivHouseholderQR ( 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
FullPivHouseholderQR()

Definition at line 93 of file FullPivHouseholderQR.h.

94  : m_qr(rows, cols),
95  m_hCoeffs((std::min)(rows,cols)),
96  m_rows_transpositions(rows),
97  m_cols_transpositions(cols),
98  m_cols_permutation(cols),
99  m_temp((std::min)(rows,cols)),
100  m_isInitialized(false),
101  m_usePrescribedThreshold(false) {}
template<typename _MatrixType>
Eigen::FullPivHouseholderQR< _MatrixType >::FullPivHouseholderQR ( 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:

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

Definition at line 115 of file FullPivHouseholderQR.h.

116  : m_qr(matrix.rows(), matrix.cols()),
117  m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
118  m_rows_transpositions(matrix.rows()),
119  m_cols_transpositions(matrix.cols()),
120  m_cols_permutation(matrix.cols()),
121  m_temp((std::min)(matrix.rows(), matrix.cols())),
122  m_isInitialized(false),
123  m_usePrescribedThreshold(false)
124  {
125  compute(matrix);
126  }
Definition: math3d.h:219
FullPivHouseholderQR & compute(const MatrixType &matrix)

Member Function Documentation

template<typename MatrixType >
MatrixType::RealScalar Eigen::FullPivHouseholderQR< 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 383 of file FullPivHouseholderQR.h.

384 {
385  using std::abs;
386  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
387  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
388  return abs(m_qr.diagonal().prod());
389 }
template<typename _MatrixType>
const PermutationType& Eigen::FullPivHouseholderQR< _MatrixType >::colsPermutation ( ) const
inline
Returns
a const reference to the column permutation matrix

Definition at line 168 of file FullPivHouseholderQR.h.

169  {
170  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
171  return m_cols_permutation;
172  }
template<typename MatrixType >
FullPivHouseholderQR< MatrixType > & Eigen::FullPivHouseholderQR< 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 FullPivHouseholderQR, FullPivHouseholderQR(const MatrixType&)

Definition at line 406 of file FullPivHouseholderQR.h.

407 {
408  using std::abs;
409  Index rows = matrix.rows();
410  Index cols = matrix.cols();
411  Index size = (std::min)(rows,cols);
412 
413  m_qr = matrix;
414  m_hCoeffs.resize(size);
415 
416  m_temp.resize(cols);
417 
418  m_precision = NumTraits<Scalar>::epsilon() * size;
419 
420  m_rows_transpositions.resize(matrix.rows());
421  m_cols_transpositions.resize(matrix.cols());
422  Index number_of_transpositions = 0;
423 
424  RealScalar biggest(0);
425 
426  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
427  m_maxpivot = RealScalar(0);
428 
429  for (Index k = 0; k < size; ++k)
430  {
431  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
432  RealScalar biggest_in_corner;
433 
434  biggest_in_corner = m_qr.bottomRightCorner(rows-k, cols-k)
435  .cwiseAbs()
436  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
437  row_of_biggest_in_corner += k;
438  col_of_biggest_in_corner += k;
439  if(k==0) biggest = biggest_in_corner;
440 
441  // if the corner is negligible, then we have less than full rank, and we can finish early
442  if(internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision))
443  {
444  m_nonzero_pivots = k;
445  for(Index i = k; i < size; i++)
446  {
447  m_rows_transpositions.coeffRef(i) = i;
448  m_cols_transpositions.coeffRef(i) = i;
449  m_hCoeffs.coeffRef(i) = Scalar(0);
450  }
451  break;
452  }
453 
454  m_rows_transpositions.coeffRef(k) = row_of_biggest_in_corner;
455  m_cols_transpositions.coeffRef(k) = col_of_biggest_in_corner;
456  if(k != row_of_biggest_in_corner) {
457  m_qr.row(k).tail(cols-k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols-k));
458  ++number_of_transpositions;
459  }
460  if(k != col_of_biggest_in_corner) {
461  m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
462  ++number_of_transpositions;
463  }
464 
465  RealScalar beta;
466  m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
467  m_qr.coeffRef(k,k) = beta;
468 
469  // remember the maximum absolute value of diagonal coefficients
470  if(abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
471 
472  m_qr.bottomRightCorner(rows-k, cols-k-1)
473  .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
474  }
475 
476  m_cols_permutation.setIdentity(cols);
477  for(Index k = 0; k < size; ++k)
478  m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
479 
480  m_det_pq = (number_of_transpositions%2) ? -1 : 1;
481  m_isInitialized = true;
482 
483  return *this;
484 }
Definition: math3d.h:219
EIGEN_STRONG_INLINE void resize(Index nbRows, Index nbCols)
Derived & applyTranspositionOnTheRight(Index i, Index j)
template<typename _MatrixType>
Index Eigen::FullPivHouseholderQR< _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 233 of file FullPivHouseholderQR.h.

234  {
235  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
236  return cols() - rank();
237  }
template<typename _MatrixType>
const HCoeffsType& Eigen::FullPivHouseholderQR< _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 297 of file FullPivHouseholderQR.h.

297 { return m_hCoeffs; }
template<typename _MatrixType>
const internal::solve_retval<FullPivHouseholderQR, typename MatrixType::IdentityReturnType> Eigen::FullPivHouseholderQR< _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 283 of file FullPivHouseholderQR.h.

284  {
285  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
286  return internal::solve_retval<FullPivHouseholderQR,typename MatrixType::IdentityReturnType>
287  (*this, MatrixType::Identity(m_qr.rows(), m_qr.cols()));
288  }
template<typename _MatrixType>
bool Eigen::FullPivHouseholderQR< _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 246 of file FullPivHouseholderQR.h.

247  {
248  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
249  return rank() == cols();
250  }
template<typename _MatrixType>
bool Eigen::FullPivHouseholderQR< _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 271 of file FullPivHouseholderQR.h.

272  {
273  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
274  return isInjective() && isSurjective();
275  }
template<typename _MatrixType>
bool Eigen::FullPivHouseholderQR< _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 259 of file FullPivHouseholderQR.h.

260  {
261  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
262  return rank() == rows();
263  }
template<typename MatrixType >
MatrixType::RealScalar Eigen::FullPivHouseholderQR< 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 392 of file FullPivHouseholderQR.h.

393 {
394  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
395  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
396  return m_qr.diagonal().cwiseAbs().array().log().sum();
397 }
template<typename MatrixType >
FullPivHouseholderQR< MatrixType >::MatrixQReturnType Eigen::FullPivHouseholderQR< MatrixType >::matrixQ ( void  ) const
inline
Returns
Expression object representing the matrix Q

Definition at line 604 of file FullPivHouseholderQR.h.

605 {
606  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
607  return MatrixQReturnType(m_qr, m_hCoeffs, m_rows_transpositions);
608 }
template<typename _MatrixType>
const MatrixType& Eigen::FullPivHouseholderQR< _MatrixType >::matrixQR ( ) const
inline
Returns
a reference to the matrix where the Householder QR decomposition is stored

Definition at line 159 of file FullPivHouseholderQR.h.

160  {
161  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
162  return m_qr;
163  }
template<typename _MatrixType>
RealScalar Eigen::FullPivHouseholderQR< _MatrixType >::maxPivot ( ) const
inline
Returns
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of U.

Definition at line 366 of file FullPivHouseholderQR.h.

366 { return m_maxpivot; }
template<typename _MatrixType>
Index Eigen::FullPivHouseholderQR< _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 357 of file FullPivHouseholderQR.h.

358  {
359  eigen_assert(m_isInitialized && "LU is not initialized.");
360  return m_nonzero_pivots;
361  }
template<typename _MatrixType>
Index Eigen::FullPivHouseholderQR< _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 216 of file FullPivHouseholderQR.h.

217  {
218  using std::abs;
219  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
220  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
221  Index result = 0;
222  for(Index i = 0; i < m_nonzero_pivots; ++i)
223  result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold);
224  return result;
225  }
template<typename _MatrixType>
const IntColVectorType& Eigen::FullPivHouseholderQR< _MatrixType >::rowsTranspositions ( ) const
inline
Returns
a const reference to the vector of indices representing the rows transpositions

Definition at line 175 of file FullPivHouseholderQR.h.

176  {
177  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
178  return m_rows_transpositions;
179  }
template<typename _MatrixType>
FullPivHouseholderQR& Eigen::FullPivHouseholderQR< _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 316 of file FullPivHouseholderQR.h.

317  {
318  m_usePrescribedThreshold = true;
319  m_prescribedThreshold = threshold;
320  return *this;
321  }
template<typename _MatrixType>
FullPivHouseholderQR& Eigen::FullPivHouseholderQR< _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 331 of file FullPivHouseholderQR.h.

332  {
333  m_usePrescribedThreshold = false;
334  return *this;
335  }
template<typename _MatrixType>
template<typename Rhs >
const internal::solve_retval<FullPivHouseholderQR, Rhs> Eigen::FullPivHouseholderQR< _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 147 of file FullPivHouseholderQR.h.

148  {
149  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
150  return internal::solve_retval<FullPivHouseholderQR, Rhs>(*this, b.derived());
151  }
template<typename _MatrixType>
RealScalar Eigen::FullPivHouseholderQR< _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 341 of file FullPivHouseholderQR.h.

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

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