44 template<
typename _MatrixType> 
class HouseholderQR
    55     typedef typename MatrixType::Scalar 
Scalar;
    96     template<
typename InputType>
   114     template<
typename InputType>
   116       : 
m_qr(matrix.derived()),
   138     template<
typename Rhs>
   169     template<
typename InputType>
   214     #ifndef EIGEN_PARSED_BY_DOXYGEN   215     template<
typename RhsType, 
typename DstType>
   217     void _solve_impl(
const RhsType &rhs, DstType &dst) 
const;
   235 template<
typename MatrixType>
   240   eigen_assert(
m_qr.rows() == 
m_qr.cols() && 
"You can't take the determinant of a non-square matrix!");
   241   return abs(
m_qr.diagonal().prod());
   244 template<
typename MatrixType>
   248   eigen_assert(
m_qr.rows() == 
m_qr.cols() && 
"You can't take the determinant of a non-square matrix!");
   249   return m_qr.diagonal().cwiseAbs().array().log().sum();
   255 template<
typename MatrixQR, 
typename HCoeffs>
   258   typedef typename MatrixQR::Scalar 
Scalar;
   259   typedef typename MatrixQR::RealScalar 
RealScalar;
   271     tempData = tempVector.data();
   274   for(
Index k = 0; k < size; ++k)
   276     Index remainingRows = rows - k;
   277     Index remainingCols = cols - k - 1;
   280     mat.col(k).tail(remainingRows).makeHouseholderInPlace(hCoeffs.coeffRef(k), beta);
   281     mat.coeffRef(k,k) = beta;
   284     mat.bottomRightCorner(remainingRows, remainingCols)
   285         .applyHouseholderOnTheLeft(mat.col(k).tail(remainingRows-1), hCoeffs.coeffRef(k), tempData+k+1);
   290 template<
typename MatrixQR, 
typename HCoeffs,
   291   typename MatrixQRScalar = 
typename MatrixQR::Scalar,
   292   bool InnerStrideIsOne = (MatrixQR::InnerStrideAtCompileTime == 1 && HCoeffs::InnerStrideAtCompileTime == 1)>
   297       typename MatrixQR::Scalar* tempData = 0)
   299     typedef typename MatrixQR::Scalar 
Scalar;
   311       tempData = tempVector.data();
   317     for (k = 0; k < size; k += blockSize)
   320       Index tcols = cols - k - bs;              
   321       Index brows = rows-k;                     
   331       BlockType A11_21 = mat.block(k,k,brows,bs);
   338         BlockType A21_22 = mat.block(k,k+bs,brows,tcols);
   347 #ifndef EIGEN_PARSED_BY_DOXYGEN   348 template<
typename _MatrixType>
   349 template<
typename RhsType, 
typename DstType>
   355   typename RhsType::PlainObject c(rhs);
   363   m_qr.topLeftCorner(rank, rank)
   364       .template triangularView<Upper>()
   365       .solveInPlace(c.topRows(rank));
   367   dst.topRows(rank) = c.topRows(rank);
   368   dst.bottomRows(
cols()-rank).setZero();
   378 template<
typename MatrixType>
   400 template<
typename Derived>
 HouseholderQR & compute(const EigenBase< InputType > &matrix)
HouseholderSequence< VectorsType, CoeffsType > householderSequence(const VectorsType &v, const CoeffsType &h)
Convenience function for constructing a Householder sequence. 
static void run(MatrixQR &mat, HCoeffs &hCoeffs, Index maxBlockSize=32, typename MatrixQR::Scalar *tempData=0)
MatrixType::StorageIndex StorageIndex
HouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation. 
void apply_block_householder_on_the_left(MatrixType &mat, const VectorsType &vectors, const CoeffsType &hCoeffs, bool forward)
const MatrixType & matrixQR() const
void householder_qr_inplace_unblocked(MatrixQR &mat, HCoeffs &hCoeffs, typename MatrixQR::Scalar *tempData=0)
Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime,(MatrixType::Flags &RowMajorBit) ? RowMajor :ColMajor, MaxRowsAtCompileTime, MaxRowsAtCompileTime > MatrixQType
static void check_template_parameters()
MatrixType::RealScalar absDeterminant() const
const Solve< HouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Sequence of Householder reflections acting on subspaces with decreasing size. 
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
MatrixType::RealScalar RealScalar
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API. 
HouseholderQR(EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix. 
MatrixType::RealScalar logAbsDeterminant() const
const HouseholderQR< PlainObject > householderQr() const
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
HouseholderSequenceType householderQ() const
internal::plain_row_type< MatrixType >::type RowVectorType
Expression of a fixed-size or dynamic-size block. 
Householder QR decomposition of a matrix. 
internal::plain_diag_type< MatrixType >::type HCoeffsType
Pseudo expression representing a solving operation. 
HouseholderQR()
Default Constructor. 
The matrix class, also used for vectors and row-vectors. 
HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::type > HouseholderSequenceType
EIGEN_DEVICE_FUNC Derived & derived()
EIGEN_DEVICE_FUNC void _solve_impl(const RhsType &rhs, DstType &dst) const
Base class for all dense matrices, vectors, and expressions. 
MatrixType::Scalar Scalar
HouseholderQR(const EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix. 
const HCoeffsType & hCoeffs() const