compbio
Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
Eigen::HouseholderQR< _MatrixType > Class Template Reference

Householder QR decomposition of a matrix. More...

#include <HouseholderQR.h>

Public Types

enum  { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime,(MatrixType::Flags &RowMajorBit) ? RowMajor :ColMajor, MaxRowsAtCompileTime, MaxRowsAtCompileTime > MatrixQType
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
typedef HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::type > HouseholderSequenceType
 

Public Member Functions

 HouseholderQR ()
 Default Constructor. More...
 
 HouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
template<typename InputType >
 HouseholderQR (const EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix. More...
 
template<typename InputType >
 HouseholderQR (EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix. More...
 
template<typename Rhs >
const Solve< HouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 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. More...
 
HouseholderSequenceType householderQ () const
 This method returns an expression of the unitary matrix Q as a sequence of Householder transformations. More...
 
const MatrixType & matrixQR () const
 
template<typename InputType >
HouseholderQRcompute (const EigenBase< InputType > &matrix)
 
MatrixType::RealScalar absDeterminant () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
Index rows () const
 
Index cols () const
 
const HCoeffsTypehCoeffs () const
 
template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 

Protected Member Functions

void computeInPlace ()
 Performs the QR factorization of the given matrix matrix. More...
 

Static Protected Member Functions

static void check_template_parameters ()
 

Protected Attributes

MatrixType m_qr
 
HCoeffsType m_hCoeffs
 
RowVectorType m_temp
 
bool m_isInitialized
 

Detailed Description

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

Householder QR decomposition of a matrix.

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

This class supports the inplace decomposition mechanism.

See also
MatrixBase::householderQr()

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);
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.
See also
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&)

§ 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 ( 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:

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.
See also
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.

§ 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.

Example:

Output:

 

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