![]() |
|
Householder QR decomposition of a matrix. More...
#include <HouseholderQR.h>
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... | |
HouseholderSequenceType | householderQ () const |
const MatrixType & | matrixQR () const |
template<typename InputType > | |
HouseholderQR & | compute (const EigenBase< InputType > &matrix) |
MatrixType::Scalar | determinant () const |
MatrixType::RealScalar | absDeterminant () const |
MatrixType::RealScalar | logAbsDeterminant () const |
MatrixType::Scalar | signDeterminant () const |
Index | rows () const |
Index | cols () const |
const HCoeffsType & | hCoeffs () const |
template<typename RhsType , typename DstType > | |
void | _solve_impl (const RhsType &rhs, DstType &dst) const |
template<bool Conjugate, typename RhsType , typename DstType > | |
void | _solve_impl_transposed (const RhsType &rhs, DstType &dst) const |
![]() | |
SolverBase () | |
~SolverBase () | |
const Solve< HouseholderQR< MatrixType_ >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
const ConstTransposeReturnType | transpose () const |
const AdjointReturnType | adjoint () const |
constexpr EIGEN_DEVICE_FUNC HouseholderQR< MatrixType_ > & | derived () |
constexpr EIGEN_DEVICE_FUNC const HouseholderQR< MatrixType_ > & | derived () const |
![]() | |
constexpr EIGEN_DEVICE_FUNC Derived & | derived () |
constexpr EIGEN_DEVICE_FUNC const Derived & | derived () const |
EIGEN_DEVICE_FUNC Derived & | const_cast_derived () const |
EIGEN_DEVICE_FUNC const Derived & | const_derived () const |
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index | size () const EIGEN_NOEXCEPT |
template<typename Dest > | |
EIGEN_DEVICE_FUNC void | evalTo (Dest &dst) const |
template<typename Dest > | |
EIGEN_DEVICE_FUNC void | addTo (Dest &dst) const |
template<typename Dest > | |
EIGEN_DEVICE_FUNC void | subTo (Dest &dst) const |
template<typename Dest > | |
EIGEN_DEVICE_FUNC void | applyThisOnTheRight (Dest &dst) const |
template<typename Dest > | |
EIGEN_DEVICE_FUNC void | applyThisOnTheLeft (Dest &dst) const |
template<typename Device > | |
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DeviceWrapper< Derived, Device > | device (Device &device) |
template<typename Device > | |
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DeviceWrapper< const Derived, Device > | device (Device &device) const |
Protected Member Functions | |
void | computeInPlace () |
![]() | |
void | _check_solve_assertion (const Rhs &b) const |
Protected Attributes | |
MatrixType | m_qr |
HCoeffsType | m_hCoeffs |
RowVectorType | m_temp |
bool | m_isInitialized |
Friends | |
class | SolverBase< HouseholderQR > |
Householder QR decomposition of a matrix.
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.
typedef SolverBase<HouseholderQR> Eigen::HouseholderQR< MatrixType_ >::Base |
typedef internal::plain_diag_type<MatrixType>::type Eigen::HouseholderQR< MatrixType_ >::HCoeffsType |
typedef HouseholderSequence<MatrixType, internal::remove_all_t<typename HCoeffsType::ConjugateReturnType> > Eigen::HouseholderQR< MatrixType_ >::HouseholderSequenceType |
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, (MatrixType::Flags & RowMajorBit) ? RowMajor : ColMajor, MaxRowsAtCompileTime, MaxRowsAtCompileTime> Eigen::HouseholderQR< MatrixType_ >::MatrixQType |
typedef MatrixType_ Eigen::HouseholderQR< MatrixType_ >::MatrixType |
typedef internal::plain_row_type<MatrixType>::type Eigen::HouseholderQR< MatrixType_ >::RowVectorType |
|
inline |
Default Constructor.
The default constructor is useful in cases in which the user intends to perform decompositions via HouseholderQR::compute(const MatrixType&).
|
inline |
Default Constructor with memory preallocation.
Like the default constructor but with preallocation of the internal data according to the specified problem size.
|
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:
References Eigen::HouseholderQR< MatrixType_ >::compute(), and matrix().
|
inlineexplicit |
Constructs a QR factorization from a given matrix.
This overloaded constructor is provided for inplace decomposition when MatrixType
is a Eigen::Ref.
References Eigen::HouseholderQR< MatrixType_ >::computeInPlace().
void Eigen::HouseholderQR< MatrixType_ >::_solve_impl | ( | const RhsType & | rhs, |
DstType & | dst | ||
) | const |
References Eigen::SolverBase< HouseholderQR< MatrixType_ > >::adjoint(), calibrate::c, Eigen::HouseholderQR< MatrixType_ >::cols(), Eigen::HouseholderQR< MatrixType_ >::householderQ(), Eigen::HouseholderQR< MatrixType_ >::m_qr, min, and Eigen::HouseholderQR< MatrixType_ >::rows().
void Eigen::HouseholderQR< MatrixType_ >::_solve_impl_transposed | ( | const RhsType & | rhs, |
DstType & | dst | ||
) | const |
MatrixType::RealScalar Eigen::HouseholderQR< MatrixType >::absDeterminant |
References abs(), eigen_assert, Eigen::HouseholderQR< MatrixType_ >::m_isInitialized, and Eigen::HouseholderQR< MatrixType_ >::m_qr.
|
inline |
References Eigen::HouseholderQR< MatrixType_ >::m_qr.
Referenced by gdb.printers._MatrixEntryIterator::__next__(), Eigen::HouseholderQR< MatrixType_ >::_solve_impl(), Eigen::HouseholderQR< MatrixType_ >::_solve_impl_transposed(), gdb.printers.EigenMatrixPrinter::children(), gdb.printers.EigenSparseMatrixPrinter::children(), Eigen::HouseholderQR< MatrixType_ >::computeInPlace(), Eigen::internal::householder_qr_inplace_unblocked(), Eigen::internal::householder_qr_inplace_blocked< MatrixQR, HCoeffs, MatrixQRScalar, InnerStrideIsOne >::run(), gdb.printers.EigenMatrixPrinter::to_string(), and gdb.printers.EigenSparseMatrixPrinter::to_string().
|
inline |
|
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.
References Eigen::HouseholderQR< MatrixType_ >::cols(), Eigen::HouseholderQR< MatrixType_ >::m_hCoeffs, Eigen::HouseholderQR< MatrixType_ >::m_isInitialized, Eigen::HouseholderQR< MatrixType_ >::m_qr, Eigen::HouseholderQR< MatrixType_ >::m_temp, min, Eigen::HouseholderQR< MatrixType_ >::rows(), Eigen::internal::householder_qr_inplace_blocked< MatrixQR, HCoeffs, MatrixQRScalar, InnerStrideIsOne >::run(), and Eigen::EigenBase< Derived >::size().
Referenced by Eigen::HouseholderQR< MatrixType_ >::compute(), and Eigen::HouseholderQR< MatrixType_ >::HouseholderQR().
MatrixType::Scalar Eigen::HouseholderQR< MatrixType >::determinant |
References eigen_assert, Eigen::HouseholderQR< MatrixType_ >::m_hCoeffs, Eigen::HouseholderQR< MatrixType_ >::m_isInitialized, Eigen::HouseholderQR< MatrixType_ >::m_qr, and Eigen::run().
|
inline |
Q
.For advanced uses only.
References Eigen::HouseholderQR< MatrixType_ >::m_hCoeffs.
Referenced by Eigen::internal::householder_qr_inplace_unblocked(), Eigen::internal::householder_qr_inplace_update(), Eigen::internal::householder_determinant< HCoeffs, Scalar, IsComplex >::run(), Eigen::internal::householder_determinant< HCoeffs, Scalar, false >::run(), and Eigen::internal::householder_qr_inplace_blocked< MatrixQR, HCoeffs, MatrixQRScalar, InnerStrideIsOne >::run().
|
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:
References eigen_assert, Eigen::HouseholderQR< MatrixType_ >::m_hCoeffs, Eigen::HouseholderQR< MatrixType_ >::m_isInitialized, and Eigen::HouseholderQR< MatrixType_ >::m_qr.
Referenced by Eigen::HouseholderQR< MatrixType_ >::_solve_impl(), Eigen::HouseholderQR< MatrixType_ >::_solve_impl_transposed(), Eigen::createRandomPIMatrixOfRank(), Eigen::RealQZ< MatrixType_ >::hessenbergTriangular(), qr(), randomMatrixWithRealEivals(), randomMatrixWithImagEivals< MatrixType, 0 >::run(), randomMatrixWithImagEivals< MatrixType, 1 >::run(), Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveNumericalDiffOneStep(), and Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveOneStep().
MatrixType::RealScalar Eigen::HouseholderQR< MatrixType >::logAbsDeterminant |
References eigen_assert, Eigen::HouseholderQR< MatrixType_ >::m_isInitialized, and Eigen::HouseholderQR< MatrixType_ >::m_qr.
|
inline |
References eigen_assert, Eigen::HouseholderQR< MatrixType_ >::m_isInitialized, and Eigen::HouseholderQR< MatrixType_ >::m_qr.
Referenced by Eigen::RealQZ< MatrixType_ >::hessenbergTriangular(), qr(), Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveNumericalDiffOneStep(), and Eigen::HybridNonLinearSolver< FunctorType, Scalar >::solveOneStep().
|
inline |
References Eigen::HouseholderQR< MatrixType_ >::m_qr.
Referenced by gdb.printers._MatrixEntryIterator::__next__(), Eigen::HouseholderQR< MatrixType_ >::_solve_impl(), Eigen::HouseholderQR< MatrixType_ >::_solve_impl_transposed(), gdb.printers.EigenMatrixPrinter::children(), gdb.printers.EigenSparseMatrixPrinter::children(), Eigen::HouseholderQR< MatrixType_ >::computeInPlace(), Eigen::internal::householder_qr_inplace_unblocked(), Eigen::internal::householder_qr_inplace_update(), Eigen::internal::householder_qr_inplace_blocked< MatrixQR, HCoeffs, MatrixQRScalar, InnerStrideIsOne >::run(), gdb.printers.EigenMatrixPrinter::to_string(), and gdb.printers.EigenSparseMatrixPrinter::to_string().
MatrixType::Scalar Eigen::HouseholderQR< MatrixType >::signDeterminant |
References eigen_assert, Eigen::HouseholderQR< MatrixType_ >::m_hCoeffs, Eigen::HouseholderQR< MatrixType_ >::m_isInitialized, Eigen::HouseholderQR< MatrixType_ >::m_qr, and Eigen::run().
|
friend |
|
protected |
|
protected |
Referenced by Eigen::HouseholderQR< MatrixType_ >::absDeterminant(), Eigen::HouseholderQR< MatrixType_ >::computeInPlace(), Eigen::HouseholderQR< MatrixType_ >::determinant(), Eigen::HouseholderQR< MatrixType_ >::householderQ(), Eigen::HouseholderQR< MatrixType_ >::logAbsDeterminant(), Eigen::HouseholderQR< MatrixType_ >::matrixQR(), and Eigen::HouseholderQR< MatrixType_ >::signDeterminant().
|
protected |
Referenced by Eigen::HouseholderQR< MatrixType_ >::_solve_impl(), Eigen::HouseholderQR< MatrixType_ >::_solve_impl_transposed(), Eigen::HouseholderQR< MatrixType_ >::absDeterminant(), Eigen::HouseholderQR< MatrixType_ >::cols(), Eigen::HouseholderQR< MatrixType_ >::compute(), Eigen::HouseholderQR< MatrixType_ >::computeInPlace(), Eigen::HouseholderQR< MatrixType_ >::determinant(), Eigen::HouseholderQR< MatrixType_ >::householderQ(), Eigen::HouseholderQR< MatrixType_ >::logAbsDeterminant(), Eigen::HouseholderQR< MatrixType_ >::matrixQR(), Eigen::HouseholderQR< MatrixType_ >::rows(), and Eigen::HouseholderQR< MatrixType_ >::signDeterminant().
|
protected |
Referenced by Eigen::HouseholderQR< MatrixType_ >::computeInPlace().