Yoji KURODA / Eigen

Dependents:   Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more

Embed: (wiki syntax)

« Back to documentation index

ColPivHouseholderQR< _MatrixType > Class Template Reference

ColPivHouseholderQR< _MatrixType > Class Template Reference
[QR module]

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

#include <ColPivHouseholderQR.h>

Public Member Functions

 ColPivHouseholderQR ()
 Default Constructor.
 ColPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation.
 ColPivHouseholderQR (const MatrixType &matrix)
 Constructs a QR factorization from a given matrix.
template<typename Rhs >
const internal::solve_retval
< ColPivHouseholderQR, Rhs > 
solve (const MatrixBase< Rhs > &b) const
 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.
HouseholderSequenceType householderQ (void) const
const MatrixType & matrixQR () const
const MatrixType & matrixR () const
ColPivHouseholderQRcompute (const MatrixType &matrix)
 Performs the QR factorization of the given matrix 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
const HCoeffsType & hCoeffs () const
ColPivHouseholderQRsetThreshold (const RealScalar &threshold)
 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.
ColPivHouseholderQRsetThreshold (Default_t)
 Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold.
RealScalar threshold () const
 Returns the threshold that will be used by certain methods such as rank().
Index nonzeroPivots () const
RealScalar maxPivot () const
ComputationInfo info () const
 Reports whether the QR factorization was succesful.

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 37 of file ColPivHouseholderQR.h.


Constructor & Destructor Documentation

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.

ColPivHouseholderQR ( Index  rows,
Index  cols 
)

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 88 of file ColPivHouseholderQR.h.

ColPivHouseholderQR ( const MatrixType &  matrix )

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 110 of file ColPivHouseholderQR.h.


Member Function Documentation

MatrixType::RealScalar 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 406 of file ColPivHouseholderQR.h.

const PermutationType& colsPermutation (  ) const
Returns:
a const reference to the column permutation matrix

Definition at line 180 of file ColPivHouseholderQR.h.

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 429 of file ColPivHouseholderQR.h.

Index dimensionOfKernel (  ) const
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 238 of file ColPivHouseholderQR.h.

const HCoeffsType& hCoeffs (  ) const
Returns:
a const reference to the vector of Householder coefficients used to represent the factor Q.

For advanced uses only.

Definition at line 303 of file ColPivHouseholderQR.h.

ColPivHouseholderQR< MatrixType >::HouseholderSequenceType householderQ ( void   ) const
Returns:
the matrix Q as a sequence of householder transformations. You can extract the meaningful part only by using:
 qr.householderQ().setLength(qr.nonzeroPivots()) 

Definition at line 561 of file ColPivHouseholderQR.h.

ComputationInfo info (  ) const

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 380 of file ColPivHouseholderQR.h.

const internal::solve_retval<ColPivHouseholderQR, typename MatrixType::IdentityReturnType> inverse (  ) const
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 289 of file ColPivHouseholderQR.h.

bool isInjective (  ) const
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 251 of file ColPivHouseholderQR.h.

bool isInvertible (  ) const
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 276 of file ColPivHouseholderQR.h.

bool isSurjective (  ) const
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 264 of file ColPivHouseholderQR.h.

MatrixType::RealScalar 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 415 of file ColPivHouseholderQR.h.

const MatrixType& matrixQR (  ) const
Returns:
a reference to the matrix where the Householder QR decomposition is stored

Definition at line 156 of file ColPivHouseholderQR.h.

const MatrixType& matrixR (  ) const
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 171 of file ColPivHouseholderQR.h.

RealScalar maxPivot (  ) const
Returns:
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of R.

Definition at line 372 of file ColPivHouseholderQR.h.

Index nonzeroPivots (  ) const
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 363 of file ColPivHouseholderQR.h.

Index rank (  ) const
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 221 of file ColPivHouseholderQR.h.

ColPivHouseholderQR& setThreshold ( Default_t   )

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 337 of file ColPivHouseholderQR.h.

ColPivHouseholderQR& setThreshold ( const RealScalar &  threshold )

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 322 of file ColPivHouseholderQR.h.

const internal::solve_retval<ColPivHouseholderQR, Rhs> solve ( const MatrixBase< Rhs > &  b ) const

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 142 of file ColPivHouseholderQR.h.

RealScalar threshold (  ) const

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

See the documentation of setThreshold(const RealScalar&).

Definition at line 347 of file ColPivHouseholderQR.h.