template<typename _MatrixType>
Eigen::HouseholderQR class

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.

Base classes

template<typename Derived>
class SolverBase
A base class for matrix decomposition and solvers.

Public types

enum (anonymous) { MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
using Base = SolverBase<HouseholderQR>
using HCoeffsType = internal::plain_diag_type<MatrixType>::type
using HouseholderSequenceType = HouseholderSequence<MatrixType, typename internal::remove_all<typename HCoeffsType::ConjugateReturnType>::type>
using MatrixQType = Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,(MatrixType::Flags&RowMajorBit) ? RowMajor :ColMajor, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
using MatrixType = _MatrixType
using RowVectorType = internal::plain_row_type<MatrixType>::type

Constructors, destructors, conversion operators

HouseholderQR()
Default Constructor.
HouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
template<typename InputType>
HouseholderQR(const EigenBase<InputType>& matrix) explicit
Constructs a QR factorization from a given matrix.
template<typename InputType>
HouseholderQR(EigenBase<InputType>& matrix) explicit
Constructs a QR factorization from a given matrix.

Public functions

auto absDeterminant() const -> MatrixType::RealScalar
auto cols() const -> Index
template<typename InputType>
auto compute(const EigenBase<InputType>& matrix) -> HouseholderQR&
auto hCoeffs() const -> const HCoeffsType&
auto householderQ() const -> HouseholderSequenceType
auto logAbsDeterminant() const -> MatrixType::RealScalar
auto matrixQR() const -> const MatrixType&
auto rows() const -> Index
template<typename Rhs>
auto solve(const MatrixBase<Rhs>& b) const -> const Solve<HouseholderQR, Rhs>

Protected static functions

static void check_template_parameters()

Protected functions

void computeInPlace()

Protected variables

HCoeffsType m_hCoeffs
bool m_isInitialized
MatrixType m_qr
RowVectorType m_temp

Enum documentation

template<typename _MatrixType>
enum Eigen::HouseholderQR<_MatrixType>::(anonymous)

Enumerators
MaxRowsAtCompileTime
MaxColsAtCompileTime

Typedef documentation

template<typename _MatrixType>
typedef SolverBase<HouseholderQR> Eigen::HouseholderQR<_MatrixType>::Base

template<typename _MatrixType>
typedef internal::plain_diag_type<MatrixType>::type Eigen::HouseholderQR<_MatrixType>::HCoeffsType

template<typename _MatrixType>
typedef HouseholderSequence<MatrixType, typename internal::remove_all<typename HCoeffsType::ConjugateReturnType>::type> Eigen::HouseholderQR<_MatrixType>::HouseholderSequenceType

template<typename _MatrixType>
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,(MatrixType::Flags&RowMajorBit) ? RowMajor :ColMajor, MaxRowsAtCompileTime, MaxRowsAtCompileTime> Eigen::HouseholderQR<_MatrixType>::MatrixQType

template<typename _MatrixType>
typedef _MatrixType Eigen::HouseholderQR<_MatrixType>::MatrixType

template<typename _MatrixType>
typedef internal::plain_row_type<MatrixType>::type Eigen::HouseholderQR<_MatrixType>::RowVectorType

Function documentation

template<typename _MatrixType>
Eigen::HouseholderQR<_MatrixType>::HouseholderQR()

Default Constructor.

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

template<typename _MatrixType>
Eigen::HouseholderQR<_MatrixType>::HouseholderQR(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.

template<typename _MatrixType> template<typename InputType>
Eigen::HouseholderQR<_MatrixType>::HouseholderQR(const EigenBase<InputType>& matrix) explicit

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

template<typename _MatrixType> template<typename InputType>
Eigen::HouseholderQR<_MatrixType>::HouseholderQR(EigenBase<InputType>& matrix) explicit

Constructs a QR factorization from a given matrix.

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

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.

template<typename _MatrixType>
Index Eigen::HouseholderQR<_MatrixType>::cols() const

template<typename _MatrixType> template<typename InputType>
HouseholderQR& Eigen::HouseholderQR<_MatrixType>::compute(const EigenBase<InputType>& matrix)

template<typename _MatrixType>
const HCoeffsType& Eigen::HouseholderQR<_MatrixType>::hCoeffs() const

Returns a const reference to the vector of Householder coefficients used to represent the factor Q.

For advanced uses only.

template<typename _MatrixType>
HouseholderSequenceType Eigen::HouseholderQR<_MatrixType>::householderQ() const

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;
A.setRandom();
HouseholderQR<MatrixXf> qr(A);
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";

Output:

The complete unitary matrix Q is:
  -0.676   0.0793    0.713  -0.0788   -0.147
  -0.221   -0.322    -0.37   -0.366   -0.759
  -0.353   -0.345   -0.214    0.841  -0.0518
   0.582   -0.462    0.555    0.176   -0.329
  -0.174   -0.747 -0.00907   -0.348    0.539

The thin matrix Q is:
  -0.676   0.0793    0.713
  -0.221   -0.322    -0.37
  -0.353   -0.345   -0.214
   0.582   -0.462    0.555
  -0.174   -0.747 -0.00907

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.

template<typename _MatrixType>
const MatrixType& Eigen::HouseholderQR<_MatrixType>::matrixQR() const

Returns a reference to the matrix where the Householder QR decomposition is stored in a LAPACK-compatible way.

template<typename _MatrixType>
Index Eigen::HouseholderQR<_MatrixType>::rows() const

template<typename _MatrixType> template<typename Rhs>
const Solve<HouseholderQR, Rhs> Eigen::HouseholderQR<_MatrixType>::solve(const MatrixBase<Rhs>& b) const

Parameters
b the right-hand-side of the equation to solve.
Returns a solution.

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.

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;
Matrix3f x;
x = m.householderQr().solve(y);
assert(y.isApprox(m*x));
cout << "Here is a solution x to the equation mx=y:" << endl << x << endl;

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

template<typename _MatrixType>
static void Eigen::HouseholderQR<_MatrixType>::check_template_parameters() protected

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.

Variable documentation

template<typename _MatrixType>
HCoeffsType Eigen::HouseholderQR<_MatrixType>::m_hCoeffs protected

template<typename _MatrixType>
bool Eigen::HouseholderQR<_MatrixType>::m_isInitialized protected

template<typename _MatrixType>
MatrixType Eigen::HouseholderQR<_MatrixType>::m_qr protected

template<typename _MatrixType>
RowVectorType Eigen::HouseholderQR<_MatrixType>::m_temp protected