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

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

#include <FullPivHouseholderQR.h>

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

Public Types

enum  { MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef MatrixType_ MatrixType
 
typedef SolverBase< FullPivHouseholderQRBase
 
typedef PermutationIndex_ PermutationIndex
 
typedef internal::FullPivHouseholderQRMatrixQReturnType< MatrixType, PermutationIndexMatrixQReturnType
 
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 PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndexPermutationType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
typedef internal::plain_col_type< MatrixType >::type ColVectorType
 
typedef MatrixType::PlainObject PlainObject
 
- Public Types inherited from Eigen::SolverBase< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >
enum  
 
typedef EigenBase< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > > Base
 
typedef internal::traits< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >::Scalar Scalar
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const FullPivHouseholderQR< 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

 FullPivHouseholderQR ()
 Default Constructor. More...
 
 FullPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation. 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...
 
MatrixQReturnType matrixQ (void) const
 
const MatrixTypematrixQR () const
 
template<typename InputType >
FullPivHouseholderQRcompute (const EigenBase< InputType > &matrix)
 
const PermutationTypecolsPermutation () const
 
const IntDiagSizeVectorTyperowsTranspositions () 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< FullPivHouseholderQRinverse () const
 
Index rows () const
 
Index cols () const
 
const HCoeffsTypehCoeffs () const
 
FullPivHouseholderQRsetThreshold (const RealScalar &threshold)
 
FullPivHouseholderQRsetThreshold (Default_t)
 
RealScalar threshold () const
 
Index nonzeroPivots () const
 
RealScalar maxPivot () const
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<bool Conjugate, typename RhsType , typename DstType >
void _solve_impl_transposed (const RhsType &rhs, DstType &dst) const
 
template<typename InputType >
FullPivHouseholderQR< MatrixType, PermutationIndex > & compute (const EigenBase< InputType > &matrix)
 
- Public Member Functions inherited from Eigen::SolverBase< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >
 SolverBase ()
 
 ~SolverBase ()
 
const Solve< FullPivHouseholderQR< MatrixType_, PermutationIndex_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const ConstTransposeReturnType transpose () const
 
const AdjointReturnType adjoint () const
 
constexpr EIGEN_DEVICE_FUNC FullPivHouseholderQR< MatrixType_, PermutationIndex_ > & derived ()
 
constexpr EIGEN_DEVICE_FUNC const FullPivHouseholderQR< 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< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >
void _check_solve_assertion (const Rhs &b) const
 

Protected Attributes

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

Friends

class SolverBase< FullPivHouseholderQR >
 

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

Member Typedef Documentation

◆ Base

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

◆ ColVectorType

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

◆ HCoeffsType

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

◆ IntDiagSizeVectorType

◆ MatrixQReturnType

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

◆ MatrixType

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

◆ PermutationIndex

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

◆ PermutationType

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

◆ PlainObject

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

◆ RowVectorType

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

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , typename PermutationIndex_ >
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
71  {
72  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
73  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
74  };
@ MaxColsAtCompileTime
Definition: FullPivHouseholderQR.h:73
@ MaxRowsAtCompileTime
Definition: FullPivHouseholderQR.h:72

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

91  : m_qr(),
92  m_hCoeffs(),
96  m_temp(),
97  m_isInitialized(false),
98  m_usePrescribedThreshold(false) {}
IntDiagSizeVectorType m_rows_transpositions
Definition: FullPivHouseholderQR.h:425
RowVectorType m_temp
Definition: FullPivHouseholderQR.h:428
IntDiagSizeVectorType m_cols_transpositions
Definition: FullPivHouseholderQR.h:426
bool m_usePrescribedThreshold
Definition: FullPivHouseholderQR.h:429
HCoeffsType m_hCoeffs
Definition: FullPivHouseholderQR.h:424
MatrixType m_qr
Definition: FullPivHouseholderQR.h:423
PermutationType m_cols_permutation
Definition: FullPivHouseholderQR.h:427
bool m_isInitialized
Definition: FullPivHouseholderQR.h:429

◆ 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()
107  : m_qr(rows, cols),
112  m_temp(cols),
113  m_isInitialized(false),
114  m_usePrescribedThreshold(false) {}
Index rows() const
Definition: FullPivHouseholderQR.h:336
Index cols() const
Definition: FullPivHouseholderQR.h:337
#define min(a, b)
Definition: datatypes.h:22

◆ 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)
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()
130  : m_qr(matrix.rows(), matrix.cols()),
131  m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
132  m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
133  m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())),
134  m_cols_permutation(matrix.cols()),
135  m_temp(matrix.cols()),
136  m_isInitialized(false),
137  m_usePrescribedThreshold(false) {
138  compute(matrix.derived());
139  }
FullPivHouseholderQR & compute(const EigenBase< InputType > &matrix)

References Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::compute(), and 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&)
150  : m_qr(matrix.derived()),
151  m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
152  m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
153  m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())),
154  m_cols_permutation(matrix.cols()),
155  m_temp(matrix.cols()),
156  m_isInitialized(false),
157  m_usePrescribedThreshold(false) {
158  computeInPlace();
159  }
void computeInPlace()
Definition: FullPivHouseholderQR.h:485

References Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::computeInPlace().

Member Function Documentation

◆ _solve_impl()

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename RhsType , typename DstType >
void Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const
564  {
565  const Index l_rank = rank();
566 
567  // FIXME introduce nonzeroPivots() and use it here. and more generally,
568  // make the same improvements in this dec as in FullPivLU.
569  if (l_rank == 0) {
570  dst.setZero();
571  return;
572  }
573 
574  typename RhsType::PlainObject c(rhs);
575 
576  Matrix<typename RhsType::Scalar, 1, RhsType::ColsAtCompileTime> temp(rhs.cols());
577  for (Index k = 0; k < l_rank; ++k) {
578  Index remainingSize = rows() - k;
579  c.row(k).swap(c.row(m_rows_transpositions.coeff(k)));
580  c.bottomRightCorner(remainingSize, rhs.cols())
581  .applyHouseholderOnTheLeft(m_qr.col(k).tail(remainingSize - 1), m_hCoeffs.coeff(k), &temp.coeffRef(0));
582  }
583 
584  m_qr.topLeftCorner(l_rank, l_rank).template triangularView<Upper>().solveInPlace(c.topRows(l_rank));
585 
586  for (Index i = 0; i < l_rank; ++i) dst.row(m_cols_permutation.indices().coeff(i)) = c.row(i);
587  for (Index i = l_rank; i < cols(); ++i) dst.row(m_cols_permutation.indices().coeff(i)).setZero();
588 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Index rank() const
Definition: FullPivHouseholderQR.h:271
const IndicesType & indices() const
Definition: PermutationMatrix.h:334
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
Definition: PlainObjectBase.h:198
res setZero()
char char char int int * k
Definition: level2_impl.h:374
int c
Definition: calibrate.py:100
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43

References calibrate::c, Eigen::PlainObjectBase< Derived >::coeff(), Eigen::Matrix< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ >::coeffRef(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols(), i, Eigen::PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, StorageIndex_ >::indices(), k, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_cols_permutation, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_rows_transpositions, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::rank(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows(), and setZero().

◆ _solve_impl_transposed()

template<typename MatrixType_ , typename PermutationIndex_ >
template<bool Conjugate, typename RhsType , typename DstType >
void Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const
593  {
594  const Index l_rank = rank();
595 
596  if (l_rank == 0) {
597  dst.setZero();
598  return;
599  }
600 
601  typename RhsType::PlainObject c(m_cols_permutation.transpose() * rhs);
602 
603  m_qr.topLeftCorner(l_rank, l_rank)
604  .template triangularView<Upper>()
605  .transpose()
606  .template conjugateIf<Conjugate>()
607  .solveInPlace(c.topRows(l_rank));
608 
609  dst.topRows(l_rank) = c.topRows(l_rank);
610  dst.bottomRows(rows() - l_rank).setZero();
611 
612  Matrix<Scalar, 1, DstType::ColsAtCompileTime> temp(dst.cols());
613  const Index size = (std::min)(rows(), cols());
614  for (Index k = size - 1; k >= 0; --k) {
615  Index remainingSize = rows() - k;
616 
617  dst.bottomRightCorner(remainingSize, dst.cols())
618  .applyHouseholderOnTheLeft(m_qr.col(k).tail(remainingSize - 1).template conjugateIf<!Conjugate>(),
619  m_hCoeffs.template conjugateIf<Conjugate>().coeff(k), &temp.coeffRef(0));
620 
621  dst.row(k).swap(dst.row(m_rows_transpositions.coeff(k)));
622  }
623 }
InverseReturnType transpose() const
Definition: PermutationMatrix.h:177
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:64

References calibrate::c, Eigen::PlainObjectBase< Derived >::coeff(), Eigen::Matrix< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ >::coeffRef(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols(), k, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_cols_permutation, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_rows_transpositions, min, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::rank(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows(), Eigen::EigenBase< Derived >::size(), and Eigen::PermutationBase< Derived >::transpose().

◆ 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()
446  {
447  using std::abs;
448  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
449  eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
450  return isInjective() ? abs(m_qr.diagonal().prod()) : RealScalar(0);
451 }
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: FullPivHouseholderQR.h:298

References abs(), eigen_assert, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::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
196  {
197  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
198  return m_cols_permutation;
199  }

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

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

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

◆ computeInPlace()

template<typename MatrixType , typename PermutationIndex >
void Eigen::FullPivHouseholderQR< MatrixType, PermutationIndex >::computeInPlace
protected
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  m_hCoeffs.resize(size);
493 
494  m_temp.resize(cols);
495 
497 
500  Index number_of_transpositions = 0;
501 
502  RealScalar biggest(0);
503 
504  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
505  m_maxpivot = RealScalar(0);
506 
507  for (Index k = 0; k < size; ++k) {
508  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
509  typedef internal::scalar_score_coeff_op<Scalar> Scoring;
510  typedef typename Scoring::result_type Score;
511 
512  Score score = m_qr.bottomRightCorner(rows - k, cols - k)
513  .unaryExpr(Scoring())
514  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
515  row_of_biggest_in_corner += k;
516  col_of_biggest_in_corner += k;
517  RealScalar biggest_in_corner =
518  internal::abs_knowing_score<Scalar>()(m_qr(row_of_biggest_in_corner, col_of_biggest_in_corner), score);
519  if (k == 0) biggest = biggest_in_corner;
520 
521  // if the corner is negligible, then we have less than full rank, and we can finish early
522  if (internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision)) {
524  for (Index i = k; i < size; i++) {
525  m_rows_transpositions.coeffRef(i) = internal::convert_index<PermutationIndex>(i);
526  m_cols_transpositions.coeffRef(i) = internal::convert_index<PermutationIndex>(i);
527  m_hCoeffs.coeffRef(i) = Scalar(0);
528  }
529  break;
530  }
531 
532  m_rows_transpositions.coeffRef(k) = internal::convert_index<PermutationIndex>(row_of_biggest_in_corner);
533  m_cols_transpositions.coeffRef(k) = internal::convert_index<PermutationIndex>(col_of_biggest_in_corner);
534  if (k != row_of_biggest_in_corner) {
535  m_qr.row(k).tail(cols - k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols - k));
536  ++number_of_transpositions;
537  }
538  if (k != col_of_biggest_in_corner) {
539  m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
540  ++number_of_transpositions;
541  }
542 
544  m_qr.col(k).tail(rows - k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
545  m_qr.coeffRef(k, k) = beta;
546 
547  // remember the maximum absolute value of diagonal coefficients
548  if (abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
549 
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 
556 
557  m_det_p = (number_of_transpositions % 2) ? -1 : 1;
558  m_isInitialized = true;
559 }
Index m_det_p
Definition: FullPivHouseholderQR.h:433
RealScalar m_precision
Definition: FullPivHouseholderQR.h:432
Index m_nonzero_pivots
Definition: FullPivHouseholderQR.h:431
RealScalar m_maxpivot
Definition: FullPivHouseholderQR.h:430
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:217
Derived & applyTranspositionOnTheRight(Index i, Index j)
Definition: PermutationMatrix.h:162
void setIdentity()
Definition: PermutationMatrix.h:122
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
internal::traits< FullPivHouseholderQR< MatrixType_, PermutationIndex_ > >::Scalar Scalar
Definition: SolverBase.h:75
Scalar beta
Definition: level2_cplx_impl.h:36
EIGEN_DEVICE_FUNC bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Definition: MathFunctions.h:1916
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43

References abs(), Eigen::PermutationBase< Derived >::applyTranspositionOnTheRight(), beta, Eigen::PlainObjectBase< Derived >::coeff(), Eigen::Matrix< Scalar_, Rows_, Cols_, Options_, MaxRows_, MaxCols_ >::coeffRef(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::cols(), eigen_assert, oomph::SarahBL::epsilon, i, Eigen::internal::isMuchSmallerThan(), k, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_cols_permutation, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_cols_transpositions, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_det_p, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_hCoeffs, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_maxpivot, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_nonzero_pivots, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_precision, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_rows_transpositions, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_temp, min, Eigen::PlainObjectBase< Derived >::resize(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::rows(), Eigen::PermutationBase< Derived >::setIdentity(), and Eigen::EigenBase< Derived >::size().

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

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

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

◆ 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&).
286  {
287  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
288  return cols() - rank();
289  }

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

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

343 { return m_hCoeffs; }

References Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::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.
331  {
332  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
333  return Inverse<FullPivHouseholderQR>(*this);
334  }

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

◆ 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&).
298  {
299  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
300  return rank() == cols();
301  }

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

Referenced by Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::absDeterminant(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::determinant(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInvertible(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::logAbsDeterminant(), and Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::signDeterminant().

◆ 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&).
321  {
322  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
323  return isInjective() && isSurjective();
324  }
bool isSurjective() const
Definition: FullPivHouseholderQR.h:310

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

◆ 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&).
310  {
311  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
312  return rank() == rows();
313  }

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

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

◆ 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()
454  {
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 isInjective() ? m_qr.diagonal().cwiseAbs().array().log().sum() : -NumTraits<RealScalar>::infinity();
458 }

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

◆ 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
704  {
705  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
707 }
internal::FullPivHouseholderQRMatrixQReturnType< MatrixType, PermutationIndex > MatrixQReturnType
Definition: FullPivHouseholderQR.h:75

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

◆ 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
187  {
188  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
189  return m_qr;
190  }

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

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

References Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::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()
400  {
401  eigen_assert(m_isInitialized && "LU is not initialized.");
402  return m_nonzero_pivots;
403  }

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

◆ 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&).
271  {
272  using std::abs;
273  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
274  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
275  Index result = 0;
276  for (Index i = 0; i < m_nonzero_pivots; ++i) result += (abs(m_qr.coeff(i, i)) > premultiplied_threshold);
277  return result;
278  }
RealScalar threshold() const
Definition: FullPivHouseholderQR.h:385

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

Referenced by Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::_solve_impl_transposed(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::dimensionOfKernel(), Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::isInjective(), and Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::isSurjective().

◆ 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
202  {
203  eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
204  return m_rows_transpositions;
205  }

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

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

362  {
365  return *this;
366  }
RealScalar m_prescribedThreshold
Definition: FullPivHouseholderQR.h:430

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

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

See the documentation of setThreshold(const RealScalar&).

376  {
377  m_usePrescribedThreshold = false;
378  return *this;
379  }

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

◆ signDeterminant()

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

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

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

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

References eigen_assert, oomph::SarahBL::epsilon, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_prescribedThreshold, Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_qr, and Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold.

Referenced by Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::rank(), and Eigen::FullPivHouseholderQR< MatrixType_, PermutationIndex_ >::setThreshold().

Friends And Related Function Documentation

◆ SolverBase< FullPivHouseholderQR >

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

Member Data Documentation

◆ m_cols_permutation

◆ m_cols_transpositions

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

◆ m_det_p

◆ m_hCoeffs

◆ m_isInitialized

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

◆ m_maxpivot

◆ m_nonzero_pivots

◆ m_precision

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

◆ m_prescribedThreshold

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

◆ m_qr

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

◆ m_rows_transpositions

◆ m_temp

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

◆ m_usePrescribedThreshold

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

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