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

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

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

Public Types

enum  {
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef SolverBase< FullPivHouseholderQRBase
 
typedef internal::plain_col_type< MatrixType >::type ColVectorType
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef Matrix< PermutationIndex, 1, internal::min_size_prefer_dynamic(ColsAtCompileTime, RowsAtCompileTime), RowMajor, 1, internal::min_size_prefer_fixed(MaxColsAtCompileTime, MaxRowsAtCompileTime)> IntDiagSizeVectorType
 
typedef internal::FullPivHouseholderQRMatrixQReturnType< MatrixType, PermutationIndexMatrixQReturnType
 
typedef MatrixType_ MatrixType
 
typedef PermutationIndex_ PermutationIndex
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndexPermutationType
 
typedef MatrixType::PlainObject PlainObject
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
- Public Types inherited from Eigen::SolverBase< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >
enum  
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
typedef EigenBase< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > > Base
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const FullPivHouseholderQR< MatrixType_, PermutationIndex_ > > ConstTransposeReturnType
 
typedef internal::traits< FullPivHouseholderQR< 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
 
Index cols () const
 
const PermutationTypecolsPermutation () const
 
template<typename InputType >
FullPivHouseholderQRcompute (const EigenBase< InputType > &matrix)
 
template<typename InputType >
FullPivHouseholderQR< MatrixType, PermutationIndex > & compute (const EigenBase< InputType > &matrix)
 
MatrixType::Scalar determinant () const
 
Index dimensionOfKernel () const
 
 FullPivHouseholderQR ()
 Default Constructor. More...
 
template<typename InputType >
 FullPivHouseholderQR (const EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix. More...
 
template<typename InputType >
 FullPivHouseholderQR (EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix. More...
 
 FullPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
const HCoeffsTypehCoeffs () const
 
const Inverse< FullPivHouseholderQRinverse () const
 
bool isInjective () const
 
bool isInvertible () const
 
bool isSurjective () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
MatrixQReturnType matrixQ (void) const
 
const MatrixTypematrixQR () const
 
RealScalar maxPivot () const
 
Index nonzeroPivots () const
 
Index rank () const
 
Index rows () const
 
const IntDiagSizeVectorTyperowsTranspositions () const
 
FullPivHouseholderQRsetThreshold (const RealScalar &threshold)
 
FullPivHouseholderQRsetThreshold (Default_t)
 
template<typename Rhs >
const Solve< FullPivHouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
- Public Member Functions inherited from Eigen::SolverBase< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >
const AdjointReturnType adjoint () const
 
FullPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived ()
 
const FullPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived () const
 
const Solve< FullPivHouseholderQR< 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< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >
void _check_solve_assertion (const Rhs &b) const
 

Protected Attributes

PermutationType m_cols_permutation
 
IntDiagSizeVectorType m_cols_transpositions
 
Index m_det_p
 
HCoeffsType m_hCoeffs
 
bool m_isInitialized
 
RealScalar m_maxpivot
 
Index m_nonzero_pivots
 
RealScalar m_precision
 
RealScalar m_prescribedThreshold
 
MatrixType m_qr
 
IntDiagSizeVectorType m_rows_transpositions
 
RowVectorType m_temp
 
bool m_usePrescribedThreshold
 

Detailed Description

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

Householder rank-revealing QR decomposition of a matrix with full 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, P', Q and R such that

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

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

This decomposition performs a very prudent full pivoting in order to be rank-revealing and achieve optimal numerical stability. The trade-off is that it is slower than HouseholderQR and ColPivHouseholderQR.

This class supports the inplace decomposition mechanism.

See also
MatrixBase::fullPivHouseholderQr()

Definition at line 62 of file FullPivHouseholderQR.h.

Member Typedef Documentation

◆ Base

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

Definition at line 68 of file FullPivHouseholderQR.h.

◆ ColVectorType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::plain_col_type<MatrixType>::type Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::ColVectorType

Definition at line 84 of file FullPivHouseholderQR.h.

◆ HCoeffsType

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

Definition at line 78 of file FullPivHouseholderQR.h.

◆ IntDiagSizeVectorType

◆ MatrixQReturnType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::FullPivHouseholderQRMatrixQReturnType<MatrixType, PermutationIndex> Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::MatrixQReturnType

Definition at line 77 of file FullPivHouseholderQR.h.

◆ MatrixType

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

Definition at line 67 of file FullPivHouseholderQR.h.

◆ PermutationIndex

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

Definition at line 70 of file FullPivHouseholderQR.h.

◆ PermutationType

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

Definition at line 82 of file FullPivHouseholderQR.h.

◆ PlainObject

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

Definition at line 85 of file FullPivHouseholderQR.h.

◆ RowVectorType

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

Definition at line 83 of file FullPivHouseholderQR.h.

Member Enumeration Documentation

◆ anonymous enum

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

Definition at line 73 of file FullPivHouseholderQR.h.

73  {
74  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
75  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
76  };

Constructor & Destructor Documentation

◆ FullPivHouseholderQR() [1/4]

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

Default Constructor.

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

Definition at line 92 of file FullPivHouseholderQR.h.

◆ FullPivHouseholderQR() [2/4]

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

Definition at line 108 of file FullPivHouseholderQR.h.

109  : m_qr(rows, cols),
114  m_temp(cols),
115  m_isInitialized(false),
116  m_usePrescribedThreshold(false) {}
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684

◆ FullPivHouseholderQR() [3/4]

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

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

Definition at line 131 of file FullPivHouseholderQR.h.

132  : m_qr(matrix.rows(), matrix.cols()),
133  m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
134  m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
135  m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())),
136  m_cols_permutation(matrix.cols()),
137  m_temp(matrix.cols()),
138  m_isInitialized(false),
140  {
141  compute(matrix.derived());
142  }
FullPivHouseholderQR & compute(const EigenBase< InputType > &matrix)

◆ FullPivHouseholderQR() [4/4]

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

Definition at line 151 of file FullPivHouseholderQR.h.

152  : m_qr(matrix.derived()),
153  m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
154  m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
155  m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())),
156  m_cols_permutation(matrix.cols()),
157  m_temp(matrix.cols()),
158  m_isInitialized(false),
160  {
161  computeInPlace();
162  }

Member Function Documentation

◆ absDeterminant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::RealScalar Eigen::FullPivHouseholderQR< 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 452 of file FullPivHouseholderQR.h.

453 {
454  using std::abs;
455  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
456  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
457  return abs(m_qr.diagonal().prod());
458 }
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::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols ( void  ) const
inline

Definition at line 337 of file FullPivHouseholderQR.h.

337 { return m_qr.cols(); }

◆ colsPermutation()

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

Definition at line 201 of file FullPivHouseholderQR.h.

202  {
203  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
204  return m_cols_permutation;
205  }

◆ compute() [1/2]

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

◆ compute() [2/2]

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

Definition at line 476 of file FullPivHouseholderQR.h.

477 {
478  m_qr = matrix.derived();
479  computeInPlace();
480  return *this;
481 }

◆ computeInPlace()

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

Definition at line 484 of file FullPivHouseholderQR.h.

485 {
486  eigen_assert(m_qr.cols() <= NumTraits<PermutationIndex>::highest());
487  using std::abs;
488  Index rows = m_qr.rows();
489  Index cols = m_qr.cols();
490  Index size = (std::min)(rows,cols);
491 
492 
493  m_hCoeffs.resize(size);
494 
495  m_temp.resize(cols);
496 
497  m_precision = NumTraits<Scalar>::epsilon() * RealScalar(size);
498 
501  Index number_of_transpositions = 0;
502 
503  RealScalar biggest(0);
504 
505  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
506  m_maxpivot = RealScalar(0);
507 
508  for (Index k = 0; k < size; ++k)
509  {
510  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
511  typedef internal::scalar_score_coeff_op<Scalar> Scoring;
512  typedef typename Scoring::result_type Score;
513 
514  Score score = m_qr.bottomRightCorner(rows-k, cols-k)
515  .unaryExpr(Scoring())
516  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
517  row_of_biggest_in_corner += k;
518  col_of_biggest_in_corner += k;
519  RealScalar biggest_in_corner = internal::abs_knowing_score<Scalar>()(m_qr(row_of_biggest_in_corner, col_of_biggest_in_corner), score);
520  if(k==0) biggest = biggest_in_corner;
521 
522  // if the corner is negligible, then we have less than full rank, and we can finish early
523  if(internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision))
524  {
525  m_nonzero_pivots = k;
526  for(Index i = k; i < size; i++)
527  {
528  m_rows_transpositions.coeffRef(i) = internal::convert_index<PermutationIndex>(i);
529  m_cols_transpositions.coeffRef(i) = internal::convert_index<PermutationIndex>(i);
530  m_hCoeffs.coeffRef(i) = Scalar(0);
531  }
532  break;
533  }
534 
535  m_rows_transpositions.coeffRef(k) = internal::convert_index<PermutationIndex>(row_of_biggest_in_corner);
536  m_cols_transpositions.coeffRef(k) = internal::convert_index<PermutationIndex>(col_of_biggest_in_corner);
537  if(k != row_of_biggest_in_corner) {
538  m_qr.row(k).tail(cols-k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols-k));
539  ++number_of_transpositions;
540  }
541  if(k != col_of_biggest_in_corner) {
542  m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
543  ++number_of_transpositions;
544  }
545 
546  RealScalar beta;
547  m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
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  m_qr.bottomRightCorner(rows-k, cols-k-1)
554  .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
555  }
556 
558  for(Index k = 0; k < size; ++k)
560 
561  m_det_p = (number_of_transpositions%2) ? -1 : 1;
562  m_isInitialized = true;
563 }
Derived & applyTranspositionOnTheRight(Index i, Index j)
constexpr const Scalar & coeff(Index rowId, Index colId) const
constexpr void resize(Index rows, Index cols)
constexpr Scalar & coeffRef(Index rowId, Index colId)
internal::traits< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >::Scalar Scalar
Definition: SolverBase.h:75
bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69

◆ determinant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::Scalar Eigen::FullPivHouseholderQR< 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 442 of file FullPivHouseholderQR.h.

443 {
444  eigen_assert(m_isInitialized && "HouseholderQR is not initialized.");
445  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
446  Scalar detQ;
447  internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
448  return m_qr.diagonal().prod() * detQ * Scalar(m_det_p);
449 }

◆ dimensionOfKernel()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::FullPivHouseholderQR< 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 281 of file FullPivHouseholderQR.h.

282  {
283  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
284  return cols() - rank();
285  }

◆ hCoeffs()

template<typename MatrixType_ , typename PermutationIndex_ >
const HCoeffsType& Eigen::FullPivHouseholderQR< 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 343 of file FullPivHouseholderQR.h.

343 { return m_hCoeffs; }

◆ inverse()

template<typename MatrixType_ , typename PermutationIndex_ >
const Inverse<FullPivHouseholderQR> Eigen::FullPivHouseholderQR< 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 330 of file FullPivHouseholderQR.h.

331  {
332  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
333  return Inverse<FullPivHouseholderQR>(*this);
334  }

◆ isInjective()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::FullPivHouseholderQR< 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 294 of file FullPivHouseholderQR.h.

295  {
296  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
297  return rank() == cols();
298  }

◆ isInvertible()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::FullPivHouseholderQR< 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 319 of file FullPivHouseholderQR.h.

320  {
321  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
322  return isInjective() && isSurjective();
323  }

◆ isSurjective()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::FullPivHouseholderQR< 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 307 of file FullPivHouseholderQR.h.

308  {
309  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
310  return rank() == rows();
311  }

◆ logAbsDeterminant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::RealScalar Eigen::FullPivHouseholderQR< 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 461 of file FullPivHouseholderQR.h.

462 {
463  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
464  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
465  return m_qr.diagonal().cwiseAbs().array().log().sum();
466 }

◆ matrixQ()

template<typename MatrixType , typename PermutationIndex >
FullPivHouseholderQR< MatrixType, PermutationIndex >::MatrixQReturnType Eigen::FullPivHouseholderQR< MatrixType, PermutationIndex >::matrixQ ( void  ) const
inline
Returns
Expression object representing the matrix Q

Definition at line 718 of file FullPivHouseholderQR.h.

719 {
720  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
722 }
internal::FullPivHouseholderQRMatrixQReturnType< MatrixType, PermutationIndex > MatrixQReturnType

◆ matrixQR()

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

Definition at line 191 of file FullPivHouseholderQR.h.

192  {
193  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
194  return m_qr;
195  }

◆ maxPivot()

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

Definition at line 412 of file FullPivHouseholderQR.h.

412 { return m_maxpivot; }

◆ nonzeroPivots()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::FullPivHouseholderQR< 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 403 of file FullPivHouseholderQR.h.

404  {
405  eigen_assert(m_isInitialized && "LU is not initialized.");
406  return m_nonzero_pivots;
407  }

◆ rank()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::FullPivHouseholderQR< 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 264 of file FullPivHouseholderQR.h.

265  {
266  using std::abs;
267  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
268  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
269  Index result = 0;
270  for(Index i = 0; i < m_nonzero_pivots; ++i)
271  result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold);
272  return result;
273  }

◆ rows()

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

Definition at line 336 of file FullPivHouseholderQR.h.

336 { return m_qr.rows(); }

◆ rowsTranspositions()

template<typename MatrixType_ , typename PermutationIndex_ >
const IntDiagSizeVectorType& Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::rowsTranspositions ( ) const
inline
Returns
a const reference to the vector of indices representing the rows transpositions

Definition at line 208 of file FullPivHouseholderQR.h.

209  {
210  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
211  return m_rows_transpositions;
212  }

◆ setThreshold() [1/2]

template<typename MatrixType_ , typename PermutationIndex_ >
FullPivHouseholderQR& Eigen::FullPivHouseholderQR< 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 362 of file FullPivHouseholderQR.h.

363  {
366  return *this;
367  }

◆ setThreshold() [2/2]

template<typename MatrixType_ , typename PermutationIndex_ >
FullPivHouseholderQR& Eigen::FullPivHouseholderQR< 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 377 of file FullPivHouseholderQR.h.

378  {
379  m_usePrescribedThreshold = false;
380  return *this;
381  }

◆ solve()

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

Parameters
bthe right-hand-side of the equation to solve.
Returns
the exact or least-square solution if the rank is greater or equal to the number of columns of A, and an arbitrary solution otherwise.

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.fullPivHouseholderQr().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::FullPivHouseholderQR< 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 387 of file FullPivHouseholderQR.h.

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

Member Data Documentation

◆ m_cols_permutation

template<typename MatrixType_ , typename PermutationIndex_ >
PermutationType Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_cols_permutation
protected

Definition at line 432 of file FullPivHouseholderQR.h.

◆ m_cols_transpositions

template<typename MatrixType_ , typename PermutationIndex_ >
IntDiagSizeVectorType Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_cols_transpositions
protected

Definition at line 431 of file FullPivHouseholderQR.h.

◆ m_det_p

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

Definition at line 438 of file FullPivHouseholderQR.h.

◆ m_hCoeffs

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

Definition at line 429 of file FullPivHouseholderQR.h.

◆ m_isInitialized

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

Definition at line 434 of file FullPivHouseholderQR.h.

◆ m_maxpivot

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

Definition at line 435 of file FullPivHouseholderQR.h.

◆ m_nonzero_pivots

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

Definition at line 436 of file FullPivHouseholderQR.h.

◆ m_precision

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_precision
protected

Definition at line 437 of file FullPivHouseholderQR.h.

◆ m_prescribedThreshold

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

Definition at line 435 of file FullPivHouseholderQR.h.

◆ m_qr

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

Definition at line 428 of file FullPivHouseholderQR.h.

◆ m_rows_transpositions

template<typename MatrixType_ , typename PermutationIndex_ >
IntDiagSizeVectorType Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_rows_transpositions
protected

Definition at line 430 of file FullPivHouseholderQR.h.

◆ m_temp

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

Definition at line 433 of file FullPivHouseholderQR.h.

◆ m_usePrescribedThreshold

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

Definition at line 434 of file FullPivHouseholderQR.h.


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