Eigen  3.4.90 (git rev 9589cc4e7fd8e4538bedef80dd36c7738977a8be)
 
Loading...
Searching...
No Matches
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ > Class Template Reference

#include <Eigen/src/QR/ColPivHouseholderQR.h>

Detailed Description

template<typename MatrixType_, typename PermutationIndex_>
class Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >

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

Template Parameters
MatrixType_the 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.

This class supports the inplace decomposition mechanism.

See also
MatrixBase::colPivHouseholderQr()
+ Inheritance diagram for Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >:

Public Member Functions

MatrixType::RealScalar absDeterminant () const
 
 ColPivHouseholderQR ()
 Default Constructor.
 
template<typename InputType>
 ColPivHouseholderQR (const EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix.
 
template<typename InputType>
 ColPivHouseholderQR (EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix.
 
 ColPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation.
 
const PermutationTypecolsPermutation () const
 
template<typename InputType>
ColPivHouseholderQR< MatrixType, PermutationIndex > & compute (const EigenBase< InputType > &matrix)
 
MatrixType::Scalar determinant () const
 
Index dimensionOfKernel () const
 
const HCoeffsType & hCoeffs () const
 
HouseholderSequenceType householderQ () const
 
ComputationInfo info () const
 Reports whether the QR factorization was successful.
 
const Inverse< ColPivHouseholderQRinverse () const
 
bool isInjective () const
 
bool isInvertible () const
 
bool isSurjective () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
const MatrixType & matrixQR () const
 
const MatrixType & matrixR () const
 
RealScalar maxPivot () const
 
Index nonzeroPivots () const
 
Index rank () const
 
ColPivHouseholderQRsetThreshold (const RealScalar &threshold)
 
ColPivHouseholderQRsetThreshold (Default_t)
 
MatrixType::Scalar signDeterminant () const
 
template<typename Rhs>
const Solve< ColPivHouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
- Public Member Functions inherited from Eigen::SolverBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
const AdjointReturnType adjoint () const
 
constexpr ColPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived ()
 
constexpr const ColPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived () const
 
const Solve< ColPivHouseholderQR< MatrixType_, PermutationIndex_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
- Public Member Functions inherited from Eigen::EigenBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
constexpr ColPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived ()
 
constexpr const ColPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 

Additional Inherited Members

- Public Types inherited from Eigen::EigenBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
typedef Eigen::Index Index
 The interface type of indices.
 

Constructor & Destructor Documentation

◆ ColPivHouseholderQR() [1/4]

template<typename MatrixType_, typename PermutationIndex_>
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR ( )
inline

Default Constructor.

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

◆ ColPivHouseholderQR() [2/4]

template<typename MatrixType_, typename PermutationIndex_>
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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()

◆ ColPivHouseholderQR() [3/4]

template<typename MatrixType_, typename PermutationIndex_>
template<typename InputType>
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR ( const EigenBase< InputType > & matrix)
inlineexplicit

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);
ColPivHouseholderQR()
Default Constructor.
Definition ColPivHouseholderQR.h:95
See also
compute()

◆ ColPivHouseholderQR() [4/4]

template<typename MatrixType_, typename PermutationIndex_>
template<typename InputType>
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR ( EigenBase< InputType > & matrix)
inlineexplicit

Constructs a QR factorization from a given matrix.

This overloaded constructor is provided for inplace decomposition when MatrixType is a Eigen::Ref.

See also
ColPivHouseholderQR(const EigenBase&)

Member Function Documentation

◆ absDeterminant()

template<typename MatrixType, typename PermutationIndex>
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::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
determinant(), logAbsDeterminant(), MatrixBase::determinant()

◆ colsPermutation()

template<typename MatrixType_, typename PermutationIndex_>
const PermutationType & Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::colsPermutation ( ) const
inline
Returns
a const reference to the column permutation matrix

◆ compute()

template<typename MatrixType_, typename PermutationIndex_>
template<typename InputType>
ColPivHouseholderQR< MatrixType, PermutationIndex > & Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::compute ( const EigenBase< InputType > & 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&)

◆ determinant()

template<typename MatrixType, typename PermutationIndex>
MatrixType::Scalar Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::determinant ( ) const
Returns
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
absDeterminant(), logAbsDeterminant(), MatrixBase::determinant()

◆ dimensionOfKernel()

template<typename MatrixType_, typename PermutationIndex_>
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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&).

◆ hCoeffs()

template<typename MatrixType_, typename PermutationIndex_>
const HCoeffsType & Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::hCoeffs ( ) const
inline
Returns
a const reference to the vector of Householder coefficients used to represent the factor Q.

For advanced uses only.

◆ householderQ()

template<typename MatrixType, typename PermutationIndex>
ColPivHouseholderQR< MatrixType, PermutationIndex >::HouseholderSequenceType Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::householderQ ( ) 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())

◆ info()

template<typename MatrixType_, typename PermutationIndex_>
ComputationInfo Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::info ( ) const
inline

Reports whether the QR factorization was successful.

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

◆ inverse()

template<typename MatrixType_, typename PermutationIndex_>
const Inverse< ColPivHouseholderQR > Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::inverse ( ) 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.

◆ isInjective()

template<typename MatrixType_, typename PermutationIndex_>
bool Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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&).

◆ isInvertible()

template<typename MatrixType_, typename PermutationIndex_>
bool Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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&).

◆ isSurjective()

template<typename MatrixType_, typename PermutationIndex_>
bool Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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&).

◆ logAbsDeterminant()

template<typename MatrixType, typename PermutationIndex>
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::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
determinant(), absDeterminant(), MatrixBase::determinant()

◆ matrixQR()

template<typename MatrixType_, typename PermutationIndex_>
const MatrixType & Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixQR ( ) const
inline
Returns
a reference to the matrix where the Householder QR decomposition is stored

◆ matrixR()

template<typename MatrixType_, typename PermutationIndex_>
const MatrixType & Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixR ( ) 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>()
const MatrixType & matrixR() const
Definition ColPivHouseholderQR.h:183
For rank-deficient matrices, use
matrixR().topLeftCorner(rank(), rank()).template triangularView<Upper>()
Index rank() const
Definition ColPivHouseholderQR.h:261

◆ maxPivot()

template<typename MatrixType_, typename PermutationIndex_>
RealScalar Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::maxPivot ( ) const
inline
Returns
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of R.

◆ nonzeroPivots()

template<typename MatrixType_, typename PermutationIndex_>
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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()

◆ rank()

template<typename MatrixType_, typename PermutationIndex_>
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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&).

◆ setThreshold() [1/2]

template<typename MatrixType_, typename PermutationIndex_>
ColPivHouseholderQR & Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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)

◆ setThreshold() [2/2]

template<typename MatrixType_, typename PermutationIndex_>
ColPivHouseholderQR & Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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&).

◆ signDeterminant()

template<typename MatrixType, typename PermutationIndex>
MatrixType::Scalar Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::signDeterminant ( ) const
Returns
the sign 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
determinant(), absDeterminant(), logAbsDeterminant(), MatrixBase::determinant()

◆ solve()

template<typename MatrixType_, typename PermutationIndex_>
template<typename Rhs>
const Solve< ColPivHouseholderQR, Rhs > Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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.

This method just tries to find as good a solution as possible. If you want to check whether a solution exists or if it is accurate, just call this function to get a result and then compute the error of this result, or use MatrixBase::isApprox() directly, for instance like this:

bool a_solution_exists = (A*result).isApprox(b, precision);

This method avoids dividing by zero, so that the non-existence of a solution doesn't by itself mean that you'll get inf or nan values.

If there exists more than one solution, this method will arbitrarily choose one.

Example:

Matrix3f m = Matrix3f::Random();
Matrix3f y = Matrix3f::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the matrix y:" << endl << y << endl;
x = m.colPivHouseholderQr().solve(y);
assert(y.isApprox(m* x));
cout << "Here is a solution x to the equation mx=y:" << endl << x << endl;
Matrix< float, 3, 3 > Matrix3f
3×3 matrix of type float.
Definition Matrix.h:478

Output:

Here is the matrix m:
 -0.824  -0.199   0.634
  0.924  -0.237   0.334
-0.0532   0.146  -0.779
Here is the matrix y:
 0.633  -0.23 0.0919
-0.573 -0.139  0.484
 0.982  0.542  0.256
Here is a solution x to the equation mx=y:
  -1.17  -0.131 -0.0219
   -5.2   -1.21   -3.53
  -2.16  -0.914  -0.989

◆ threshold()

template<typename MatrixType_, typename PermutationIndex_>
RealScalar Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::threshold ( ) const
inline

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

See the documentation of setThreshold(const RealScalar&).


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