Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
HouseholderQR< MatrixType > Class Template Reference

Detailed Description

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

Householder QR decomposition of a matrix.

Parameters
MatrixTypethe 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.

See Also
MatrixBase::householderQr()

Public Member Functions

MatrixType::RealScalar absDeterminant () const
 
HouseholderQRcompute (const MatrixType &matrix)
 
const HCoeffsType & hCoeffs () const
 
HouseholderSequenceType householderQ () const
 
 HouseholderQR ()
 Default Constructor. More...
 
 HouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
 HouseholderQR (const MatrixType &matrix)
 Constructs a QR factorization from a given matrix. More...
 
MatrixType::RealScalar logAbsDeterminant () const
 
const MatrixType & matrixQR () const
 
template<typename Rhs >
const internal::solve_retval
< HouseholderQR, Rhs > 
solve (const MatrixBase< Rhs > &b) const
 

Constructor & Destructor Documentation

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 ( 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 ( const MatrixType &  matrix)
inline

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);
*
See Also
compute()

References HouseholderQR< MatrixType >::compute().

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

Referenced by HouseholderQR< MatrixType >::HouseholderQR().

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

For advanced uses only.

HouseholderSequenceType householderQ ( void  ) 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:

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

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()
const MatrixType& matrixQR ( ) const
inline
Returns
a reference to the matrix where the Householder QR decomposition is stored in a LAPACK-compatible way.
const internal::solve_retval<HouseholderQR, Rhs> 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.
Note
The case where b is a matrix is not yet implemented. Also, this code is space inefficient.

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

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

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