Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
Eigen::FullPivHouseholderQR< MatrixType_ > Class Template Reference

Detailed Description

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

Householder rank-revealing QR decomposition of a matrix with full 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, P', Q and R such that

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

by using Householder transformations. Here, P and P' are permutation matrices, 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.

This class supports the inplace decomposition mechanism.

See also
MatrixBase::fullPivHouseholderQr()
+ Inheritance diagram for Eigen::FullPivHouseholderQR< MatrixType_ >:

Public Member Functions

MatrixType::RealScalar absDeterminant () const
 
const PermutationTypecolsPermutation () const
 
template<typename InputType >
FullPivHouseholderQR< MatrixType > & compute (const EigenBase< InputType > &matrix)
 
Index dimensionOfKernel () const
 
 FullPivHouseholderQR ()
 Default Constructor. More...
 
template<typename InputType >
 FullPivHouseholderQR (const EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix. More...
 
template<typename InputType >
 FullPivHouseholderQR (EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix. More...
 
 FullPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
const HCoeffsType & hCoeffs () const
 
const Inverse< FullPivHouseholderQRinverse () const
 
bool isInjective () const
 
bool isInvertible () const
 
bool isSurjective () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
MatrixQReturnType matrixQ (void) const
 
const MatrixType & matrixQR () const
 
RealScalar maxPivot () const
 
Index nonzeroPivots () const
 
Index rank () const
 
const IntDiagSizeVectorTyperowsTranspositions () const
 
FullPivHouseholderQRsetThreshold (const RealScalar &threshold)
 
FullPivHouseholderQRsetThreshold (Default_t)
 
template<typename Rhs >
const Solve< FullPivHouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
- Public Member Functions inherited from Eigen::SolverBase< FullPivHouseholderQR< MatrixType_ > >
const AdjointReturnType adjoint () const
 
FullPivHouseholderQR< MatrixType_ > & derived ()
 
const FullPivHouseholderQR< MatrixType_ > & derived () const
 
const Solve< FullPivHouseholderQR< MatrixType_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & derived ()
 
const Derived & 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< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 

Constructor & Destructor Documentation

◆ FullPivHouseholderQR() [1/4]

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&).

◆ FullPivHouseholderQR() [2/4]

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()

◆ FullPivHouseholderQR() [3/4]

template<typename MatrixType_ >
template<typename InputType >
Eigen::FullPivHouseholderQR< MatrixType_ >::FullPivHouseholderQR ( 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:

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

◆ FullPivHouseholderQR() [4/4]

template<typename MatrixType_ >
template<typename InputType >
Eigen::FullPivHouseholderQR< MatrixType_ >::FullPivHouseholderQR ( 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
FullPivHouseholderQR(const EigenBase&)

Member Function Documentation

◆ absDeterminant()

template<typename MatrixType >
MatrixType::RealScalar Eigen::FullPivHouseholderQR< MatrixType >::absDeterminant
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()

◆ colsPermutation()

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

◆ compute()

template<typename MatrixType_ >
template<typename InputType >
FullPivHouseholderQR<MatrixType>& Eigen::FullPivHouseholderQR< MatrixType_ >::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 FullPivHouseholderQR, FullPivHouseholderQR(const MatrixType&)

◆ dimensionOfKernel()

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&).

◆ hCoeffs()

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.

◆ inverse()

template<typename MatrixType_ >
const Inverse<FullPivHouseholderQR> Eigen::FullPivHouseholderQR< MatrixType_ >::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_ >
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&).

◆ isInvertible()

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&).

◆ isSurjective()

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&).

◆ logAbsDeterminant()

template<typename MatrixType >
MatrixType::RealScalar Eigen::FullPivHouseholderQR< MatrixType >::logAbsDeterminant
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()

◆ matrixQ()

template<typename MatrixType >
FullPivHouseholderQR< MatrixType >::MatrixQReturnType Eigen::FullPivHouseholderQR< MatrixType >::matrixQ ( void  ) const
inline
Returns
Expression object representing the matrix Q

◆ matrixQR()

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

◆ maxPivot()

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.

◆ nonzeroPivots()

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()

◆ rank()

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&).

◆ rowsTranspositions()

template<typename MatrixType_ >
const IntDiagSizeVectorType& Eigen::FullPivHouseholderQR< MatrixType_ >::rowsTranspositions ( ) const
inline
Returns
a const reference to the vector of indices representing the rows transpositions

◆ setThreshold() [1/2]

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)

◆ setThreshold() [2/2]

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&).

◆ solve()

template<typename MatrixType_ >
template<typename Rhs >
const Solve<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.

Parameters
bthe right-hand-side of the equation to solve.
Returns
the exact or least-square solution if the rank is greater or equal to the number of columns of A, and an arbitrary solution otherwise.

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:

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

Output:

Here is the matrix m:
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
Here is the matrix y:
  0.108   -0.27   0.832
-0.0452  0.0268   0.271
  0.258   0.904   0.435
Here is a solution x to the equation mx=y:
 0.609   2.68   1.67
-0.231  -1.57 0.0713
  0.51   3.51   1.05

◆ threshold()

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&).


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