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

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

#include <ColPivHouseholderQR.h>

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

Public Types

enum  { MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef MatrixType_ MatrixType
 
typedef SolverBase< ColPivHouseholderQRBase
 
typedef PermutationIndex_ PermutationIndex
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndexPermutationType
 
typedef internal::plain_row_type< MatrixType, PermutationIndex >::type IntRowVectorType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
typedef internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
 
typedef HouseholderSequence< MatrixType, internal::remove_all_t< typename HCoeffsType::ConjugateReturnType > > HouseholderSequenceType
 
typedef MatrixType::PlainObject PlainObject
 
- Public Types inherited from Eigen::SolverBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
enum  
 
typedef EigenBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > > Base
 
typedef internal::traits< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >::Scalar Scalar
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const ColPivHouseholderQR< MatrixType_, PermutationIndex_ > > ConstTransposeReturnType
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
- 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

 ColPivHouseholderQR ()
 Default Constructor. More...
 
 ColPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. 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...
 
HouseholderSequenceType householderQ () const
 
HouseholderSequenceType matrixQ () const
 
const MatrixTypematrixQR () const
 
const MatrixTypematrixR () const
 
template<typename InputType >
ColPivHouseholderQRcompute (const EigenBase< InputType > &matrix)
 
const PermutationTypecolsPermutation () const
 
MatrixType::Scalar determinant () const
 
MatrixType::RealScalar absDeterminant () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
MatrixType::Scalar signDeterminant () const
 
Index rank () const
 
Index dimensionOfKernel () const
 
bool isInjective () const
 
bool isSurjective () const
 
bool isInvertible () const
 
const Inverse< ColPivHouseholderQRinverse () const
 
Index rows () const
 
Index cols () const
 
const HCoeffsTypehCoeffs () const
 
ColPivHouseholderQRsetThreshold (const RealScalar &threshold)
 
ColPivHouseholderQRsetThreshold (Default_t)
 
RealScalar threshold () const
 
Index nonzeroPivots () const
 
RealScalar maxPivot () const
 
ComputationInfo info () const
 Reports whether the QR factorization was successful. More...
 
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
 
template<typename InputType >
ColPivHouseholderQR< MatrixType, PermutationIndex > & compute (const EigenBase< InputType > &matrix)
 
- Public Member Functions inherited from Eigen::SolverBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
 SolverBase ()
 
 ~SolverBase ()
 
const Solve< ColPivHouseholderQR< MatrixType_, PermutationIndex_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const ConstTransposeReturnType transpose () const
 
const AdjointReturnType adjoint () const
 
constexpr EIGEN_DEVICE_FUNC ColPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived ()
 
constexpr EIGEN_DEVICE_FUNC const ColPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived () const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
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 ()
 
- Protected Member Functions inherited from Eigen::SolverBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >
void _check_solve_assertion (const Rhs &b) const
 

Protected Attributes

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

Private Member Functions

void init (Index rows, Index cols)
 

Friends

class SolverBase< ColPivHouseholderQR >
 
class CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >
 

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

Member Typedef Documentation

◆ Base

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

◆ HCoeffsType

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

◆ HouseholderSequenceType

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

◆ IntRowVectorType

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

◆ MatrixType

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

◆ PermutationIndex

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

◆ PermutationType

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

◆ PlainObject

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

◆ RealRowVectorType

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

◆ RowVectorType

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

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , typename PermutationIndex_ >
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
62  {
63  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
64  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
65  };
@ MaxRowsAtCompileTime
Definition: ColPivHouseholderQR.h:63
@ MaxColsAtCompileTime
Definition: ColPivHouseholderQR.h:64

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&).

96  : m_qr(),
97  m_hCoeffs(),
100  m_temp(),
103  m_isInitialized(false),
104  m_usePrescribedThreshold(false) {}
bool m_usePrescribedThreshold
Definition: ColPivHouseholderQR.h:433
HCoeffsType m_hCoeffs
Definition: ColPivHouseholderQR.h:427
RealRowVectorType m_colNormsUpdated
Definition: ColPivHouseholderQR.h:431
PermutationType m_colsPermutation
Definition: ColPivHouseholderQR.h:428
RowVectorType m_temp
Definition: ColPivHouseholderQR.h:430
RealRowVectorType m_colNormsDirect
Definition: ColPivHouseholderQR.h:432
MatrixType m_qr
Definition: ColPivHouseholderQR.h:426
bool m_isInitialized
Definition: ColPivHouseholderQR.h:433
IntRowVectorType m_colsTranspositions
Definition: ColPivHouseholderQR.h:429

◆ 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()
112 : m_qr(rows, cols) { init(rows, cols); }
Index rows() const
Definition: ColPivHouseholderQR.h:326
Index cols() const
Definition: ColPivHouseholderQR.h:327
void init(Index rows, Index cols)
Definition: ColPivHouseholderQR.h:76

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::init(), and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows().

◆ 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)
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
See also
compute()
127  : m_qr(matrix.rows(), matrix.cols()) {
128  init(matrix.rows(), matrix.cols());
129  compute(matrix.derived());
130  }
ColPivHouseholderQR & compute(const EigenBase< InputType > &matrix)

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::compute(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::init(), and 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&)
140  : m_qr(matrix.derived()) {
141  init(matrix.rows(), matrix.cols());
142  computeInPlace();
143  }
void computeInPlace()
Definition: ColPivHouseholderQR.h:488

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::computeInPlace(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::init(), and matrix().

Member Function Documentation

◆ _solve_impl()

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename RhsType , typename DstType >
void Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const
588  {
589  const Index nonzero_pivots = nonzeroPivots();
590 
591  if (nonzero_pivots == 0) {
592  dst.setZero();
593  return;
594  }
595 
596  typename RhsType::PlainObject c(rhs);
597 
598  c.applyOnTheLeft(householderQ().setLength(nonzero_pivots).adjoint());
599 
600  m_qr.topLeftCorner(nonzero_pivots, nonzero_pivots)
601  .template triangularView<Upper>()
602  .solveInPlace(c.topRows(nonzero_pivots));
603 
604  for (Index i = 0; i < nonzero_pivots; ++i) dst.row(m_colsPermutation.indices().coeff(i)) = c.row(i);
605  for (Index i = nonzero_pivots; i < cols(); ++i) dst.row(m_colsPermutation.indices().coeff(i)).setZero();
606 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Index nonzeroPivots() const
Definition: ColPivHouseholderQR.h:390
HouseholderSequenceType householderQ() const
Definition: ColPivHouseholderQR.h:656
const IndicesType & indices() const
Definition: PermutationMatrix.h:334
const AdjointReturnType adjoint() const
Definition: SolverBase.h:136
res setZero()
int c
Definition: calibrate.py:100
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43

References Eigen::SolverBase< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >::adjoint(), calibrate::c, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::householderQ(), i, Eigen::PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, StorageIndex_ >::indices(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colsPermutation, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::nonzeroPivots(), and setZero().

◆ _solve_impl_transposed()

template<typename MatrixType_ , typename PermutationIndex_ >
template<bool Conjugate, typename RhsType , typename DstType >
void Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const
611  {
612  const Index nonzero_pivots = nonzeroPivots();
613 
614  if (nonzero_pivots == 0) {
615  dst.setZero();
616  return;
617  }
618 
619  typename RhsType::PlainObject c(m_colsPermutation.transpose() * rhs);
620 
621  m_qr.topLeftCorner(nonzero_pivots, nonzero_pivots)
622  .template triangularView<Upper>()
623  .transpose()
624  .template conjugateIf<Conjugate>()
625  .solveInPlace(c.topRows(nonzero_pivots));
626 
627  dst.topRows(nonzero_pivots) = c.topRows(nonzero_pivots);
628  dst.bottomRows(rows() - nonzero_pivots).setZero();
629 
630  dst.applyOnTheLeft(householderQ().setLength(nonzero_pivots).template conjugateIf<!Conjugate>());
631 }
InverseReturnType transpose() const
Definition: PermutationMatrix.h:177

References calibrate::c, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::householderQ(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colsPermutation, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::nonzeroPivots(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows(), and Eigen::PermutationBase< Derived >::transpose().

◆ 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()
449  {
450  using std::abs;
451  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
452  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
453  return isInjective() ? abs(m_qr.diagonal().prod()) : RealScalar(0);
454 }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
#define eigen_assert(x)
Definition: Macros.h:910
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
bool isInjective() const
Definition: ColPivHouseholderQR.h:288

References abs(), eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr.

◆ cols()

◆ colsPermutation()

template<typename MatrixType_ , typename PermutationIndex_ >
const PermutationType& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::colsPermutation ( ) const
inline

◆ 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&)
481  {
482  m_qr = matrix.derived();
483  computeInPlace();
484  return *this;
485 }

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::computeInPlace(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, and matrix().

◆ computeInPlace()

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

References abs(), Eigen::numext::abs2(), Eigen::PermutationBase< Derived >::applyTranspositionOnTheRight(), beta, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols(), eigen_assert, Eigen::numext::is_exactly_zero(), j, k, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colNormsDirect, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colNormsUpdated, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colsPermutation, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colsTranspositions, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_det_p, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_maxpivot, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_nonzero_pivots, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_temp, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows(), Eigen::PermutationBase< Derived >::setIdentity(), Eigen::EigenBase< Derived >::size(), Eigen::numext::sqrt(), and swap().

Referenced by Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR(), and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::compute().

◆ 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()
440  {
441  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
442  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
443  Scalar detQ;
445  return isInjective() ? (detQ * Scalar(m_det_p)) * m_qr.diagonal().prod() : Scalar(0);
446 }
SCALAR Scalar
Definition: bench_gemm.cpp:45
internal::traits< ColPivHouseholderQR< MatrixType_, PermutationIndex_ > >::Scalar Scalar
Definition: SolverBase.h:75
@ IsComplex
Definition: common.h:73
auto run(Kernel kernel, Args &&... args) -> decltype(kernel(args...))
Definition: gpu_test_helper.h:414

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_det_p, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, and Eigen::run().

◆ 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&).
276  {
277  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
278  return cols() - rank();
279  }
Index rank() const
Definition: ColPivHouseholderQR.h:261

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols(), eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rank().

Referenced by Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::dimensionOfKernel().

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

333 { return m_hCoeffs; }

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs.

Referenced by Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::hCoeffs(), and Eigen::LevenbergMarquardt< FunctorType_ >::minimizeOptimumStorageOneStep().

◆ 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())
656  {
657  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
658  return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
659 }
HouseholderSequence< MatrixType, internal::remove_all_t< typename HCoeffsType::ConjugateReturnType > > HouseholderSequenceType
Definition: ColPivHouseholderQR.h:72

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr.

Referenced by Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl_transposed(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixQ(), and Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::matrixQ().

◆ 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
406  {
407  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
408  return Success;
409  }
@ Success
Definition: Constants.h:440

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::Success.

◆ init()

template<typename MatrixType_ , typename PermutationIndex_ >
void Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::init ( Index  rows,
Index  cols 
)
inlineprivate
76  {
78  m_hCoeffs.resize(diag);
80  m_colsTranspositions.resize(cols);
81  m_temp.resize(cols);
82  m_colNormsUpdated.resize(cols);
83  m_colNormsDirect.resize(cols);
84  m_isInitialized = false;
86  }
void resize(Index newSize)
Definition: PermutationMatrix.h:119
const char const char const char * diag
Definition: level2_impl.h:86
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
Definition: MathFunctions.h:920

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols(), diag, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colNormsDirect, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colNormsUpdated, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colsPermutation, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_colsTranspositions, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_temp, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold, Eigen::numext::mini(), Eigen::PermutationBase< Derived >::resize(), and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows().

Referenced by Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColPivHouseholderQR().

◆ 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.
321  {
322  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
323  return Inverse<ColPivHouseholderQR>(*this);
324  }

References eigen_assert, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized.

◆ isInjective()

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

◆ 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&).
311  {
312  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
313  return isInjective() && isSurjective();
314  }
bool isSurjective() const
Definition: ColPivHouseholderQR.h:300

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isSurjective(), and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized.

Referenced by Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::isInvertible().

◆ 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&).
300  {
301  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
302  return rank() == rows();
303  }

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rank(), and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows().

Referenced by Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInvertible(), and Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::isSurjective().

◆ 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()
457  {
458  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
459  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
460  return isInjective() ? m_qr.diagonal().cwiseAbs().array().log().sum() : -NumTraits<RealScalar>::infinity();
461 }

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr.

◆ matrixQ()

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

◆ matrixQR()

template<typename MatrixType_ , typename PermutationIndex_ >
const MatrixType& Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixQR ( ) const
inline

◆ 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
Definition: ColPivHouseholderQR.h:183
For rank-deficient matrices, use
matrixR().topLeftCorner(rank(), rank()).template triangularView<Upper>()
183  {
184  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
185  return m_qr;
186  }

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr.

◆ 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.
398 { return m_maxpivot; }

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_maxpivot.

Referenced by Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::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()
390  {
391  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
392  return m_nonzero_pivots;
393  }

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_nonzero_pivots.

Referenced by Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl_transposed(), and Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::nonzeroPivots().

◆ 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&).
261  {
262  using std::abs;
263  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
264  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
265  Index result = 0;
266  for (Index i = 0; i < m_nonzero_pivots; ++i) result += (abs(m_qr.coeff(i, i)) > premultiplied_threshold);
267  return result;
268  }
RealScalar threshold() const
Definition: ColPivHouseholderQR.h:375

References abs(), eigen_assert, i, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_maxpivot, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_nonzero_pivots, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::threshold().

Referenced by Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::dimensionOfKernel(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isSurjective(), Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::rank(), and test_sparseqr_scalar().

◆ 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)

352  {
355  return *this;
356  }
RealScalar m_prescribedThreshold
Definition: ColPivHouseholderQR.h:434

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_prescribedThreshold, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold, and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::threshold().

Referenced by Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::setThreshold().

◆ 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:361

See the documentation of setThreshold(const RealScalar&).

366  {
367  m_usePrescribedThreshold = false;
368  return *this;
369  }

References Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold.

◆ signDeterminant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::Scalar Eigen::ColPivHouseholderQR< MatrixType, PermutationIndex >::signDeterminant
Returns
the sign 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(), logAbsDeterminant(), MatrixBase::determinant()
464  {
465  eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
466  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
467  Scalar detQ;
469  return isInjective() ? (detQ * Scalar(m_det_p)) * m_qr.diagonal().array().sign().prod() : Scalar(0);
470 }

References eigen_assert, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_det_p, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, and Eigen::run().

◆ threshold()

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

Friends And Related Function Documentation

◆ CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >

template<typename MatrixType_ , typename PermutationIndex_ >
friend class CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >
friend

◆ SolverBase< ColPivHouseholderQR >

template<typename MatrixType_ , typename PermutationIndex_ >
friend class SolverBase< ColPivHouseholderQR >
friend

Member Data Documentation

◆ m_colNormsDirect

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

◆ m_colNormsUpdated

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

◆ m_colsPermutation

◆ m_colsTranspositions

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

◆ m_det_p

◆ m_hCoeffs

◆ m_isInitialized

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

Referenced by Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::_check_solve_assertion(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::absDeterminant(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::colsPermutation(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::computeInPlace(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::determinant(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::dimensionOfKernel(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::householderQ(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::info(), Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::info(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::init(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::inverse(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInvertible(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::isSurjective(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::logAbsDeterminant(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixQR(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::matrixR(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::nonzeroPivots(), Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::pseudoInverse(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::rank(), Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::signDeterminant(), and Eigen::ColPivHouseholderQR< MatrixType_, PermutationIndex_ >::threshold().

◆ m_maxpivot

◆ m_nonzero_pivots

◆ m_prescribedThreshold

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

◆ m_qr

◆ m_temp

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

◆ m_usePrescribedThreshold


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