Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
Eigen::HouseholderQR< MatrixType_ > Class Template Reference

Detailed Description

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

Householder QR decomposition of a matrix.

Template Parameters
MatrixType_the type of the matrix of which we are computing the QR decomposition

This class performs a QR decomposition of a matrix A into matrices Q and R such that

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

by using Householder transformations. Here, Q a unitary matrix and R an upper triangular matrix. The result is stored in a compact way compatible with LAPACK.

Note that no pivoting is performed. This is not a rank-revealing decomposition. If you want that feature, use FullPivHouseholderQR or ColPivHouseholderQR instead.

This Householder QR decomposition is faster, but less numerically stable and less feature-full than FullPivHouseholderQR or ColPivHouseholderQR.

This class supports the inplace decomposition mechanism.

See also
MatrixBase::householderQr()
+ Inheritance diagram for Eigen::HouseholderQR< MatrixType_ >:

Public Member Functions

MatrixType::RealScalar absDeterminant () const
 
MatrixType::Scalar determinant () const
 
const HCoeffsType & hCoeffs () const
 
HouseholderSequenceType householderQ () const
 
 HouseholderQR ()
 Default Constructor.
 
template<typename InputType >
 HouseholderQR (const EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix.
 
template<typename InputType >
 HouseholderQR (EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix.
 
 HouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation.
 
MatrixType::RealScalar logAbsDeterminant () const
 
const MatrixType & matrixQR () const
 
MatrixType::Scalar signDeterminant () const
 
template<typename Rhs >
const Solve< HouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 
- Public Member Functions inherited from Eigen::SolverBase< HouseholderQR< MatrixType_ > >
const AdjointReturnType adjoint () const
 
HouseholderQR< MatrixType_ > & derived ()
 
const HouseholderQR< MatrixType_ > & derived () const
 
const Solve< HouseholderQR< MatrixType_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
- Public Member Functions inherited from Eigen::EigenBase< HouseholderQR< MatrixType_ > >
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
HouseholderQR< MatrixType_ > & derived ()
 
const HouseholderQR< MatrixType_ > & derived () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 

Protected Member Functions

void computeInPlace ()
 

Additional Inherited Members

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

Constructor & Destructor Documentation

◆ HouseholderQR() [1/4]

template<typename MatrixType_ >
Eigen::HouseholderQR< MatrixType_ >::HouseholderQR ( )
inline

Default Constructor.

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

◆ HouseholderQR() [2/4]

template<typename MatrixType_ >
Eigen::HouseholderQR< MatrixType_ >::HouseholderQR ( 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
HouseholderQR()

◆ HouseholderQR() [3/4]

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

HouseholderQR<MatrixType> qr(matrix.rows(), matrix.cols());
qr.compute(matrix);
Householder QR decomposition of a matrix.
Definition HouseholderQR.h:59
See also
compute()

◆ HouseholderQR() [4/4]

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

Member Function Documentation

◆ absDeterminant()

template<typename MatrixType >
MatrixType::RealScalar Eigen::HouseholderQR< 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. Also, do not rely on the determinant being exactly zero for testing singularity or rank-deficiency.
See also
determinant(), logAbsDeterminant(), MatrixBase::determinant()

◆ computeInPlace()

template<typename MatrixType >
void Eigen::HouseholderQR< MatrixType >::computeInPlace ( )
protected

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 HouseholderQR, HouseholderQR(const MatrixType&)

◆ determinant()

template<typename MatrixType >
MatrixType::Scalar Eigen::HouseholderQR< MatrixType >::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. Also, do not rely on the determinant being exactly zero for testing singularity or rank-deficiency.
See also
absDeterminant(), logAbsDeterminant(), MatrixBase::determinant()

◆ hCoeffs()

template<typename MatrixType_ >
const HCoeffsType & Eigen::HouseholderQR< MatrixType_ >::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_ >
HouseholderSequenceType Eigen::HouseholderQR< MatrixType_ >::householderQ ( ) const
inline

This method returns an expression of the unitary matrix Q as a sequence of Householder transformations.

The returned expression can directly be used to perform matrix products. It can also be assigned to a dense Matrix object. Here is an example showing how to recover the full or thin matrix Q, as well as how to perform matrix products using operator*:

Example:

MatrixXf A(MatrixXf::Random(5, 3)), thinQ(MatrixXf::Identity(5, 3)), Q;
Q = qr.householderQ();
thinQ = qr.householderQ() * thinQ;
std::cout << "The complete unitary matrix Q is:\n" << Q << "\n\n";
std::cout << "The thin matrix Q is:\n" << thinQ << "\n\n";
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Derived & setRandom(Index size)
Definition Random.h:147

Output:

 

◆ logAbsDeterminant()

template<typename MatrixType >
MatrixType::RealScalar Eigen::HouseholderQR< 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.
Warning
Do not rely on the determinant being exactly zero for testing singularity or rank-deficiency.
See also
determinant(), absDeterminant(), MatrixBase::determinant()

◆ matrixQR()

template<typename MatrixType_ >
const MatrixType & Eigen::HouseholderQR< MatrixType_ >::matrixQR ( ) const
inline
Returns
a reference to the matrix where the Householder QR decomposition is stored in a LAPACK-compatible way.

◆ signDeterminant()

template<typename MatrixType >
MatrixType::Scalar Eigen::HouseholderQR< MatrixType >::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.
Warning
Do not rely on the determinant being exactly zero for testing singularity or rank-deficiency.
See also
determinant(), absDeterminant(), MatrixBase::determinant()

◆ solve()

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

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:

typedef Matrix<float, 3, 3> Matrix3x3;
Matrix3x3 m = Matrix3x3::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.householderQr().solve(y);
assert(y.isApprox(m* x));
cout << "Here is a solution x to the equation mx=y:" << endl << x << endl;
const HouseholderQR< PlainObject > householderQr() const
Definition HouseholderQR.h:526

Output:

 

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