Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ > Class Template Reference

Householder rank-revealing QR decomposition of a matrix with column-pivoting. More...

+ Inheritance diagram for Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >:

Public Types

enum  {
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef SolverBase< ColPivHouseholderQRBase
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef HouseholderSequence< MatrixType, internal::remove_all_t< typename HCoeffsType::ConjugateReturnType > > HouseholderSequenceType
 
typedef internal::plain_row_type< MatrixType, PermutationIndex >::type IntRowVectorType
 
typedef MatrixType_ MatrixType
 
typedef PermutationIndex_ PermutationIndex
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndexPermutationType
 
typedef MatrixType::PlainObject PlainObject
 
typedef internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
- Public Types inherited from Eigen::SolverBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
enum  
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
typedef EigenBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > > Base
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const ColPivHouseholderQR< MatrixType_, PermutationIndex_ > > ConstTransposeReturnType
 
typedef internal::traits< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >::Scalar Scalar
 
- Public Types inherited from Eigen::EigenBase< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 
typedef internal::traits< Derived >::StorageKind StorageKind
 

Public Member Functions

MatrixType::RealScalar absDeterminant () const
 
 ColPivHouseholderQR ()
 Default Constructor. More...
 
template<typename InputType >
 ColPivHouseholderQR (const EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix. More...
 
template<typename InputType >
 ColPivHouseholderQR (EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix. More...
 
 ColPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
Index cols () const
 
const PermutationTypecolsPermutation () const
 
template<typename InputType >
ColPivHouseholderQRcompute (const EigenBase< InputType > &matrix)
 
template<typename InputType >
ColPivHouseholderQR< MatrixType, PermutationIndex > & compute (const EigenBase< InputType > &matrix)
 
MatrixType::Scalar determinant () const
 
Index dimensionOfKernel () const
 
const HCoeffsTypehCoeffs () const
 
HouseholderSequenceType householderQ () const
 
ComputationInfo info () const
 Reports whether the QR factorization was successful. More...
 
const Inverse< ColPivHouseholderQRinverse () const
 
bool isInjective () const
 
bool isInvertible () const
 
bool isSurjective () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
HouseholderSequenceType matrixQ () const
 
const MatrixTypematrixQR () const
 
const MatrixTypematrixR () const
 
RealScalar maxPivot () const
 
Index nonzeroPivots () const
 
Index rank () const
 
Index rows () const
 
ColPivHouseholderQRsetThreshold (const RealScalar &threshold)
 
ColPivHouseholderQRsetThreshold (Default_t)
 
template<typename Rhs >
const Solve< ColPivHouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
- Public Member Functions inherited from Eigen::SolverBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
const AdjointReturnType adjoint () const
 
ColPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived ()
 
const ColPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived () const
 
const Solve< ColPivHouseholderQR< MatrixType_, PermutationIndex_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
 ~SolverBase ()
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
template<typename Dest >
void addTo (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheLeft (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheRight (Dest &dst) const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & const_cast_derived () const
 
const Derived & const_derived () const
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Dest >
void evalTo (Dest &dst) const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
template<typename Dest >
void subTo (Dest &dst) const
 

Protected Member Functions

void computeInPlace ()
 
- Protected Member Functions inherited from Eigen::SolverBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
void _check_solve_assertion (const Rhs &b) const
 

Protected Attributes

RealRowVectorType m_colNormsDirect
 
RealRowVectorType m_colNormsUpdated
 
PermutationType m_colsPermutation
 
IntRowVectorType m_colsTranspositions
 
Index m_det_p
 
HCoeffsType m_hCoeffs
 
bool m_isInitialized
 
RealScalar m_maxpivot
 
Index m_nonzero_pivots
 
RealScalar m_prescribedThreshold
 
MatrixType m_qr
 
RowVectorType m_temp
 
bool m_usePrescribedThreshold
 

Private Member Functions

void init (Index rows, Index cols)
 

Detailed Description

template<typename MatrixType_, typename PermutationIndex_>
class Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >

Householder rank-revealing QR decomposition of a matrix with column-pivoting.

Template Parameters
MatrixType_the type of the matrix of which we are computing the QR decomposition

This class performs a rank-revealing QR decomposition of a matrix A into matrices P, Q and R such that

\[ \mathbf{A} \, \mathbf{P} = \mathbf{Q} \, \mathbf{R} \]

by using Householder transformations. Here, P is a permutation matrix, Q a unitary matrix and R an upper triangular matrix.

This decomposition performs column pivoting in order to be rank-revealing and improve numerical stability. It is slower than HouseholderQR, and faster than FullPivHouseholderQR.

This class supports the inplace decomposition mechanism.

See also
MatrixBase::colPivHouseholderQr()

Definition at line 53 of file ColPivHouseholderQR.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , typename PermutationIndex_ >
typedef SolverBase<ColPivHouseholderQR> Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::Base

Definition at line 59 of file ColPivHouseholderQR.h.

◆ HCoeffsType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::plain_diag_type<MatrixType>::type Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::HCoeffsType

Definition at line 68 of file ColPivHouseholderQR.h.

◆ HouseholderSequenceType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef HouseholderSequence<MatrixType,internal::remove_all_t<typename HCoeffsType::ConjugateReturnType> > Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::HouseholderSequenceType

Definition at line 73 of file ColPivHouseholderQR.h.

◆ IntRowVectorType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::plain_row_type<MatrixType, PermutationIndex>::type Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::IntRowVectorType

Definition at line 70 of file ColPivHouseholderQR.h.

◆ MatrixType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef MatrixType_ Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::MatrixType

Definition at line 58 of file ColPivHouseholderQR.h.

◆ PermutationIndex

template<typename MatrixType_ , typename PermutationIndex_ >
typedef PermutationIndex_ Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::PermutationIndex

Definition at line 61 of file ColPivHouseholderQR.h.

◆ PermutationType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndex> Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::PermutationType

Definition at line 69 of file ColPivHouseholderQR.h.

◆ PlainObject

template<typename MatrixType_ , typename PermutationIndex_ >
typedef MatrixType::PlainObject Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::PlainObject

Definition at line 74 of file ColPivHouseholderQR.h.

◆ RealRowVectorType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::plain_row_type<MatrixType, RealScalar>::type Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::RealRowVectorType

Definition at line 72 of file ColPivHouseholderQR.h.

◆ RowVectorType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::plain_row_type<MatrixType>::type Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::RowVectorType

Definition at line 71 of file ColPivHouseholderQR.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , typename PermutationIndex_ >
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 64 of file ColPivHouseholderQR.h.

64  {
65  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
66  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
67  };

Constructor & Destructor Documentation

◆ ColPivHouseholderQR() [1/4]

template<typename MatrixType_ , typename PermutationIndex_ >
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR ( )
inline

Default Constructor.

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

Definition at line 97 of file ColPivHouseholderQR.h.

◆ ColPivHouseholderQR() [2/4]

template<typename MatrixType_ , typename PermutationIndex_ >
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR ( 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
ColPivHouseholderQR()

Definition at line 114 of file ColPivHouseholderQR.h.

114 : m_qr(rows, cols) { init(rows, cols); }
void init(Index rows, Index cols)

◆ ColPivHouseholderQR() [3/4]

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename InputType >
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR ( 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:

ColPivHouseholderQR<MatrixType> qr(matrix.rows(), matrix.cols());
qr.compute(matrix);
HouseholderQR< MatrixXf > qr(A)
See also
compute()

Definition at line 129 of file ColPivHouseholderQR.h.

129  : m_qr(matrix.rows(), matrix.cols()) {
130  init(matrix.rows(), matrix.cols());
131  compute(matrix.derived());
132  }
ColPivHouseholderQR & compute(const EigenBase< InputType > &matrix)

◆ ColPivHouseholderQR() [4/4]

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename InputType >
Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR ( 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
ColPivHouseholderQR(const EigenBase&)

Definition at line 141 of file ColPivHouseholderQR.h.

141  : m_qr(matrix.derived()) {
142  init(matrix.rows(), matrix.cols());
143  computeInPlace();
144  }

Member Function Documentation

◆ absDeterminant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::absDeterminant
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
determinant(), logAbsDeterminant(), MatrixBase::determinant()

Definition at line 457 of file ColPivHouseholderQR.h.

458 {
459  using std::abs;
460  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
461  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
462  return abs(m_qr.diagonal().prod());
463 }
const AbsReturnType abs() const
#define eigen_assert(x)
Definition: Macros.h:902
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)

◆ cols()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols ( void  ) const
inline

Definition at line 328 of file ColPivHouseholderQR.h.

328 { return m_qr.cols(); }

◆ colsPermutation()

template<typename MatrixType_ , typename PermutationIndex_ >
const PermutationType& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::colsPermutation ( ) const
inline
Returns
a const reference to the column permutation matrix

Definition at line 199 of file ColPivHouseholderQR.h.

200  {
201  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
202  return m_colsPermutation;
203  }

◆ compute() [1/2]

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename InputType >
ColPivHouseholderQR& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::compute ( const EigenBase< InputType > &  matrix)

◆ compute() [2/2]

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename InputType >
ColPivHouseholderQR<MatrixType, PermutationIndex>& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::compute ( const EigenBase< InputType > &  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 ColPivHouseholderQR, ColPivHouseholderQR(const MatrixType&)

Definition at line 481 of file ColPivHouseholderQR.h.

482 {
483  m_qr = matrix.derived();
484  computeInPlace();
485  return *this;
486 }

◆ computeInPlace()

template<typename MatrixType , typename PermutationIndex >
void Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::computeInPlace
protected

Definition at line 489 of file ColPivHouseholderQR.h.

490 {
491 
492  eigen_assert(m_qr.cols()<=NumTraits<PermutationIndex>::highest());
493 
494  using std::abs;
495 
496  Index rows = m_qr.rows();
497  Index cols = m_qr.cols();
498  Index size = m_qr.diagonalSize();
499 
500  m_hCoeffs.resize(size);
501 
502  m_temp.resize(cols);
503 
504  m_colsTranspositions.resize(m_qr.cols());
505  Index number_of_transpositions = 0;
506 
507  m_colNormsUpdated.resize(cols);
508  m_colNormsDirect.resize(cols);
509  for (Index k = 0; k < cols; ++k) {
510  // colNormsDirect(k) caches the most recent directly computed norm of
511  // column k.
512  m_colNormsDirect.coeffRef(k) = m_qr.col(k).norm();
513  m_colNormsUpdated.coeffRef(k) = m_colNormsDirect.coeffRef(k);
514  }
515 
516  RealScalar threshold_helper = numext::abs2<RealScalar>(m_colNormsUpdated.maxCoeff() * NumTraits<RealScalar>::epsilon()) / RealScalar(rows);
517  RealScalar norm_downdate_threshold = numext::sqrt(NumTraits<RealScalar>::epsilon());
518 
519  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
520  m_maxpivot = RealScalar(0);
521 
522  for(Index k = 0; k < size; ++k)
523  {
524  // first, we look up in our table m_colNormsUpdated which column has the biggest norm
525  Index biggest_col_index;
526  RealScalar biggest_col_sq_norm = numext::abs2(m_colNormsUpdated.tail(cols-k).maxCoeff(&biggest_col_index));
527  biggest_col_index += k;
528 
529  // Track the number of meaningful pivots but do not stop the decomposition to make
530  // sure that the initial matrix is properly reproduced. See bug 941.
531  if(m_nonzero_pivots==size && biggest_col_sq_norm < threshold_helper * RealScalar(rows-k))
532  m_nonzero_pivots = k;
533 
534  // apply the transposition to the columns
535  m_colsTranspositions.coeffRef(k) = static_cast<PermutationIndex>(biggest_col_index);
536  if(k != biggest_col_index) {
537  m_qr.col(k).swap(m_qr.col(biggest_col_index));
538  std::swap(m_colNormsUpdated.coeffRef(k), m_colNormsUpdated.coeffRef(biggest_col_index));
539  std::swap(m_colNormsDirect.coeffRef(k), m_colNormsDirect.coeffRef(biggest_col_index));
540  ++number_of_transpositions;
541  }
542 
543  // generate the householder vector, store it below the diagonal
544  RealScalar beta;
545  m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
546 
547  // apply the householder transformation to the diagonal coefficient
548  m_qr.coeffRef(k,k) = beta;
549 
550  // remember the maximum absolute value of diagonal coefficients
551  if(abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
552 
553  // apply the householder transformation
554  m_qr.bottomRightCorner(rows-k, cols-k-1)
555  .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
556 
557  // update our table of norms of the columns
558  for (Index j = k + 1; j < cols; ++j) {
559  // The following implements the stable norm downgrade step discussed in
560  // http://www.netlib.org/lapack/lawnspdf/lawn176.pdf
561  // and used in LAPACK routines xGEQPF and xGEQP3.
562  // See lines 278-297 in http://www.netlib.org/lapack/explore-html/dc/df4/sgeqpf_8f_source.html
563  if (!numext::is_exactly_zero(m_colNormsUpdated.coeffRef(j))) {
564  RealScalar temp = abs(m_qr.coeffRef(k, j)) / m_colNormsUpdated.coeffRef(j);
565  temp = (RealScalar(1) + temp) * (RealScalar(1) - temp);
566  temp = temp < RealScalar(0) ? RealScalar(0) : temp;
567  RealScalar temp2 = temp * numext::abs2<RealScalar>(m_colNormsUpdated.coeffRef(j) /
568  m_colNormsDirect.coeffRef(j));
569  if (temp2 <= norm_downdate_threshold) {
570  // The updated norm has become too inaccurate so re-compute the column
571  // norm directly.
572  m_colNormsDirect.coeffRef(j) = m_qr.col(j).tail(rows - k - 1).norm();
573  m_colNormsUpdated.coeffRef(j) = m_colNormsDirect.coeffRef(j);
574  } else {
575  m_colNormsUpdated.coeffRef(j) *= numext::sqrt(temp);
576  }
577  }
578  }
579  }
580 
582  for(Index k = 0; k < size/*m_nonzero_pivots*/; ++k)
584 
585  m_det_p = (number_of_transpositions%2) ? -1 : 1;
586  m_isInitialized = true;
587 }
PermutationIndex_ PermutationIndex
Derived & applyTranspositionOnTheRight(Index i, Index j)
void swap(scoped_array< T > &a, scoped_array< T > &b)
Definition: Memory.h:788
bool abs2(bool x)
EIGEN_ALWAYS_INLINE float sqrt(const float &x)
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69
std::ptrdiff_t j

◆ determinant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::Scalar Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::determinant
Returns
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
absDeterminant(), logAbsDeterminant(), MatrixBase::determinant()

Definition at line 447 of file ColPivHouseholderQR.h.

448 {
449  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
450  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
451  Scalar detQ;
452  internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
453  return m_qr.diagonal().prod() * detQ * Scalar(m_det_p);
454 }
internal::traits< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >::Scalar Scalar
Definition: SolverBase.h:75

◆ dimensionOfKernel()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::dimensionOfKernel ( ) const
inline
Returns
the dimension of the kernel of the matrix of which *this is the QR decomposition.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 272 of file ColPivHouseholderQR.h.

273  {
274  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
275  return cols() - rank();
276  }

◆ hCoeffs()

template<typename MatrixType_ , typename PermutationIndex_ >
const HCoeffsType& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::hCoeffs ( ) const
inline
Returns
a const reference to the vector of Householder coefficients used to represent the factor Q.

For advanced uses only.

Definition at line 334 of file ColPivHouseholderQR.h.

334 { return m_hCoeffs; }

◆ householderQ()

template<typename MatrixType , typename PermutationIndex >
ColPivHouseholderQR< MatrixType, PermutationIndex >::HouseholderSequenceType Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::householderQ
Returns
the matrix Q as a sequence of householder transformations. You can extract the meaningful part only by using:
qr.householderQ().setLength(qr.nonzeroPivots())

Definition at line 660 of file ColPivHouseholderQR.h.

661 {
662  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
663  return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
664 }
HouseholderSequence< MatrixType, internal::remove_all_t< typename HCoeffsType::ConjugateReturnType > > HouseholderSequenceType

◆ info()

template<typename MatrixType_ , typename PermutationIndex_ >
ComputationInfo Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::info ( ) const
inline

Reports whether the QR factorization was successful.

Note
This function always returns Success. It is provided for compatibility with other factorization routines.
Returns
Success

Definition at line 411 of file ColPivHouseholderQR.h.

412  {
413  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
414  return Success;
415  }
@ Success
Definition: Constants.h:446

◆ init()

template<typename MatrixType_ , typename PermutationIndex_ >
void Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::init ( Index  rows,
Index  cols 
)
inlineprivate

Definition at line 77 of file ColPivHouseholderQR.h.

77  {
78  Index diag = numext::mini(rows, cols);
79  m_hCoeffs.resize(diag);
81  m_colsTranspositions.resize(cols);
82  m_temp.resize(cols);
83  m_colNormsUpdated.resize(cols);
84  m_colNormsDirect.resize(cols);
85  m_isInitialized = false;
87  }
void resize(Index newSize)
EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)

◆ inverse()

template<typename MatrixType_ , typename PermutationIndex_ >
const Inverse<ColPivHouseholderQR> Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::inverse ( ) const
inline
Returns
the inverse of the matrix of which *this is the QR decomposition.
Note
If this matrix is not invertible, the returned matrix has undefined coefficients. Use isInvertible() to first determine whether this matrix is invertible.

Definition at line 321 of file ColPivHouseholderQR.h.

322  {
323  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
324  return Inverse<ColPivHouseholderQR>(*this);
325  }

◆ isInjective()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective ( ) const
inline
Returns
true if the matrix of which *this is the QR decomposition represents an injective linear map, i.e. has trivial kernel; false otherwise.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 285 of file ColPivHouseholderQR.h.

286  {
287  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
288  return rank() == cols();
289  }

◆ isInvertible()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInvertible ( ) const
inline
Returns
true if the matrix of which *this is the QR decomposition is invertible.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 310 of file ColPivHouseholderQR.h.

311  {
312  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
313  return isInjective() && isSurjective();
314  }

◆ isSurjective()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isSurjective ( ) const
inline
Returns
true if the matrix of which *this is the QR decomposition represents a surjective linear map; false otherwise.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 298 of file ColPivHouseholderQR.h.

299  {
300  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
301  return rank() == rows();
302  }

◆ logAbsDeterminant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::logAbsDeterminant
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
determinant(), absDeterminant(), MatrixBase::determinant()

Definition at line 466 of file ColPivHouseholderQR.h.

467 {
468  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
469  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
470  return m_qr.diagonal().cwiseAbs().array().log().sum();
471 }

◆ matrixQ()

template<typename MatrixType_ , typename PermutationIndex_ >
HouseholderSequenceType Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixQ ( ) const
inline

Definition at line 167 of file ColPivHouseholderQR.h.

168  {
169  return householderQ();
170  }
HouseholderSequenceType householderQ() const

◆ matrixQR()

template<typename MatrixType_ , typename PermutationIndex_ >
const MatrixType& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixQR ( ) const
inline
Returns
a reference to the matrix where the Householder QR decomposition is stored

Definition at line 174 of file ColPivHouseholderQR.h.

175  {
176  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
177  return m_qr;
178  }

◆ matrixR()

template<typename MatrixType_ , typename PermutationIndex_ >
const MatrixType& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixR ( ) const
inline
Returns
a reference to the matrix where the result Householder QR is stored
Warning
The strict lower part of this matrix contains internal values. Only the upper triangular part should be referenced. To get it, use
matrixR().template triangularView<Upper>()
const MatrixType & matrixR() const
For rank-deficient matrices, use
matrixR().topLeftCorner(rank(), rank()).template triangularView<Upper>()

Definition at line 189 of file ColPivHouseholderQR.h.

190  {
191  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
192  return m_qr;
193  }

◆ maxPivot()

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::maxPivot ( ) const
inline
Returns
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of R.

Definition at line 403 of file ColPivHouseholderQR.h.

403 { return m_maxpivot; }

◆ nonzeroPivots()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::nonzeroPivots ( ) const
inline
Returns
the number of nonzero pivots in the QR decomposition. Here nonzero is meant in the exact sense, not in a fuzzy sense. So that notion isn't really intrinsically interesting, but it is still useful when implementing algorithms.
See also
rank()

Definition at line 394 of file ColPivHouseholderQR.h.

395  {
396  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
397  return m_nonzero_pivots;
398  }

◆ rank()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rank ( ) const
inline
Returns
the rank of the matrix of which *this is the QR decomposition.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 255 of file ColPivHouseholderQR.h.

256  {
257  using std::abs;
258  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
259  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
260  Index result = 0;
261  for(Index i = 0; i < m_nonzero_pivots; ++i)
262  result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold);
263  return result;
264  }

◆ rows()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows ( void  ) const
inline

Definition at line 327 of file ColPivHouseholderQR.h.

327 { return m_qr.rows(); }

◆ setThreshold() [1/2]

template<typename MatrixType_ , typename PermutationIndex_ >
ColPivHouseholderQR& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::setThreshold ( const RealScalar &  threshold)
inline

Allows to prescribe a threshold to be used by certain methods, such as rank(), who need to determine when pivots are to be considered nonzero. This is not used for the QR decomposition itself.

When it needs to get the threshold value, Eigen calls threshold(). By default, this uses a formula to automatically determine a reasonable threshold. Once you have called the present method setThreshold(const RealScalar&), your value is used instead.

Parameters
thresholdThe new value to use as the threshold.

A pivot will be considered nonzero if its absolute value is strictly greater than \( \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert \) where maxpivot is the biggest pivot.

If you want to come back to the default behavior, call setThreshold(Default_t)

Definition at line 353 of file ColPivHouseholderQR.h.

354  {
357  return *this;
358  }

◆ setThreshold() [2/2]

template<typename MatrixType_ , typename PermutationIndex_ >
ColPivHouseholderQR& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::setThreshold ( Default_t  )
inline

Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold.

You should pass the special object Eigen::Default as parameter here.

qr.setThreshold(Eigen::Default);
@ Default
Definition: Constants.h:364

See the documentation of setThreshold(const RealScalar&).

Definition at line 368 of file ColPivHouseholderQR.h.

369  {
370  m_usePrescribedThreshold = false;
371  return *this;
372  }

◆ solve()

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename Rhs >
const Solve<ColPivHouseholderQR, Rhs> Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::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.

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);
Array< int, 3, 1 > b
MatrixXcf A
bool isApprox(const Scalar &x, const Scalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_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:

cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the matrix y:" << endl << y << endl;
x = m.colPivHouseholderQr().solve(y);
assert(y.isApprox(m*x));
cout << "Here is a solution x to the equation mx=y:" << endl << x << endl;
Matrix3f m
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< float, 3, 3 > Matrix3f
3×3 matrix of type float.
Definition: Matrix.h:501
const Scalar & y

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

◆ threshold()

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::threshold ( ) const
inline

Returns the threshold that will be used by certain methods such as rank().

See the documentation of setThreshold(const RealScalar&).

Definition at line 378 of file ColPivHouseholderQR.h.

379  {
382  // this formula comes from experimenting (see "LU precision tuning" thread on the list)
383  // and turns out to be identical to Higham's formula used already in LDLt.
384  : NumTraits<Scalar>::epsilon() * RealScalar(m_qr.diagonalSize());
385  }

Member Data Documentation

◆ m_colNormsDirect

template<typename MatrixType_ , typename PermutationIndex_ >
RealRowVectorType Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colNormsDirect
protected

Definition at line 439 of file ColPivHouseholderQR.h.

◆ m_colNormsUpdated

template<typename MatrixType_ , typename PermutationIndex_ >
RealRowVectorType Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colNormsUpdated
protected

Definition at line 438 of file ColPivHouseholderQR.h.

◆ m_colsPermutation

template<typename MatrixType_ , typename PermutationIndex_ >
PermutationType Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colsPermutation
protected

Definition at line 435 of file ColPivHouseholderQR.h.

◆ m_colsTranspositions

template<typename MatrixType_ , typename PermutationIndex_ >
IntRowVectorType Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colsTranspositions
protected

Definition at line 436 of file ColPivHouseholderQR.h.

◆ m_det_p

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_det_p
protected

Definition at line 443 of file ColPivHouseholderQR.h.

◆ m_hCoeffs

template<typename MatrixType_ , typename PermutationIndex_ >
HCoeffsType Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs
protected

Definition at line 434 of file ColPivHouseholderQR.h.

◆ m_isInitialized

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized
protected

Definition at line 440 of file ColPivHouseholderQR.h.

◆ m_maxpivot

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_maxpivot
protected

Definition at line 441 of file ColPivHouseholderQR.h.

◆ m_nonzero_pivots

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_nonzero_pivots
protected

Definition at line 442 of file ColPivHouseholderQR.h.

◆ m_prescribedThreshold

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_prescribedThreshold
protected

Definition at line 441 of file ColPivHouseholderQR.h.

◆ m_qr

template<typename MatrixType_ , typename PermutationIndex_ >
MatrixType Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr
protected

Definition at line 433 of file ColPivHouseholderQR.h.

◆ m_temp

template<typename MatrixType_ , typename PermutationIndex_ >
RowVectorType Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_temp
protected

Definition at line 437 of file ColPivHouseholderQR.h.

◆ m_usePrescribedThreshold

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold
protected

Definition at line 440 of file ColPivHouseholderQR.h.


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