Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::ColPivHouseholderQR< _MatrixType > Class Template Reference

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

#include <src/eigen/Eigen/src/QR/ColPivHouseholderQR.h>

+ Inheritance diagram for Eigen::ColPivHouseholderQR< _MatrixType >:
+ Collaboration diagram for Eigen::ColPivHouseholderQR< _MatrixType >:

Public Types

enum  { RowsAtCompileTime = MatrixType::RowsAtCompileTime , ColsAtCompileTime = MatrixType::ColsAtCompileTime , MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTimePermutationType
 
typedef internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
typedef internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
 
typedef HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::type > HouseholderSequenceType
 
typedef MatrixType::PlainObject PlainObject
 

Public Member Functions

 ColPivHouseholderQR ()
 Default Constructor.
 
 ColPivHouseholderQR (Index rows, Index cols)
 Default Constructor with memory preallocation.
 
template<typename InputType >
 ColPivHouseholderQR (const EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix.
 
template<typename InputType >
 ColPivHouseholderQR (EigenBase< InputType > &matrix)
 Constructs a QR factorization from a given matrix.
 
template<typename Rhs >
const Solve< ColPivHouseholderQR, Rhs > solve (const MatrixBase< Rhs > &b) const
 
HouseholderSequenceType householderQ () const
 
HouseholderSequenceType matrixQ () const
 
const MatrixTypematrixQR () const
 
const MatrixTypematrixR () const
 
template<typename InputType >
ColPivHouseholderQRcompute (const EigenBase< InputType > &matrix)
 
const PermutationTypecolsPermutation () const
 
MatrixType::RealScalar absDeterminant () const
 
MatrixType::RealScalar logAbsDeterminant () 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 succesful.
 
template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<typename InputType >
ColPivHouseholderQR< MatrixType > & compute (const EigenBase< InputType > &matrix)
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 

Protected Member Functions

void computeInPlace ()
 

Static Protected Member Functions

static void check_template_parameters ()
 

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_pq
 

Private Types

typedef PermutationType::StorageIndex PermIndexType
 

Friends

class CompleteOrthogonalDecomposition< MatrixType >
 

Detailed Description

template<typename _MatrixType>
class Eigen::ColPivHouseholderQR< _MatrixType >

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

Template Parameters
_MatrixTypethe 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

◆ HCoeffsType

template<typename _MatrixType >
typedef internal::plain_diag_type<MatrixType>::type Eigen::ColPivHouseholderQR< _MatrixType >::HCoeffsType

◆ HouseholderSequenceType

template<typename _MatrixType >
typedef HouseholderSequence<MatrixType,typename internal::remove_all<typename HCoeffsType::ConjugateReturnType>::type> Eigen::ColPivHouseholderQR< _MatrixType >::HouseholderSequenceType

◆ IntRowVectorType

template<typename _MatrixType >
typedef internal::plain_row_type<MatrixType,Index>::type Eigen::ColPivHouseholderQR< _MatrixType >::IntRowVectorType

◆ MatrixType

template<typename _MatrixType >
typedef _MatrixType Eigen::ColPivHouseholderQR< _MatrixType >::MatrixType

◆ PermIndexType

template<typename _MatrixType >
typedef PermutationType::StorageIndex Eigen::ColPivHouseholderQR< _MatrixType >::PermIndexType
private

◆ PermutationType

template<typename _MatrixType >
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> Eigen::ColPivHouseholderQR< _MatrixType >::PermutationType

◆ PlainObject

template<typename _MatrixType >
typedef MatrixType::PlainObject Eigen::ColPivHouseholderQR< _MatrixType >::PlainObject

◆ RealRowVectorType

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

◆ RealScalar

template<typename _MatrixType >
typedef MatrixType::RealScalar Eigen::ColPivHouseholderQR< _MatrixType >::RealScalar

◆ RowVectorType

template<typename _MatrixType >
typedef internal::plain_row_type<MatrixType>::type Eigen::ColPivHouseholderQR< _MatrixType >::RowVectorType

◆ Scalar

template<typename _MatrixType >
typedef MatrixType::Scalar Eigen::ColPivHouseholderQR< _MatrixType >::Scalar

◆ StorageIndex

template<typename _MatrixType >
typedef MatrixType::StorageIndex Eigen::ColPivHouseholderQR< _MatrixType >::StorageIndex

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType >
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
53 {
54 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
55 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
56 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
57 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
58 };
@ RowsAtCompileTime
Definition ColPivHouseholderQR.h:54
@ ColsAtCompileTime
Definition ColPivHouseholderQR.h:55
@ MaxColsAtCompileTime
Definition ColPivHouseholderQR.h:57
@ MaxRowsAtCompileTime
Definition ColPivHouseholderQR.h:56

Constructor & Destructor Documentation

◆ ColPivHouseholderQR() [1/4]

template<typename _MatrixType >
Eigen::ColPivHouseholderQR< _MatrixType >::ColPivHouseholderQR ( )
inline

Default Constructor.

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

84 : m_qr(),
85 m_hCoeffs(),
88 m_temp(),
91 m_isInitialized(false),
RowVectorType m_temp
Definition ColPivHouseholderQR.h:438
HCoeffsType m_hCoeffs
Definition ColPivHouseholderQR.h:435
bool m_isInitialized
Definition ColPivHouseholderQR.h:441
IntRowVectorType m_colsTranspositions
Definition ColPivHouseholderQR.h:437
RealRowVectorType m_colNormsUpdated
Definition ColPivHouseholderQR.h:439
bool m_usePrescribedThreshold
Definition ColPivHouseholderQR.h:441
MatrixType m_qr
Definition ColPivHouseholderQR.h:434
PermutationType m_colsPermutation
Definition ColPivHouseholderQR.h:436
RealRowVectorType m_colNormsDirect
Definition ColPivHouseholderQR.h:440

◆ ColPivHouseholderQR() [2/4]

template<typename _MatrixType >
Eigen::ColPivHouseholderQR< _MatrixType >::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()
101 : m_qr(rows, cols),
102 m_hCoeffs((std::min)(rows,cols)),
105 m_temp(cols),
108 m_isInitialized(false),
Index cols() const
Definition ColPivHouseholderQR.h:328
Index rows() const
Definition ColPivHouseholderQR.h:327
PermutationType::StorageIndex PermIndexType
Definition ColPivHouseholderQR.h:73

◆ ColPivHouseholderQR() [3/4]

template<typename _MatrixType >
template<typename InputType >
Eigen::ColPivHouseholderQR< _MatrixType >::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);
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Definition ColPivHouseholderQR.h:49
See also
compute()
125 : m_qr(matrix.rows(), matrix.cols()),
126 m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
127 m_colsPermutation(PermIndexType(matrix.cols())),
128 m_colsTranspositions(matrix.cols()),
129 m_temp(matrix.cols()),
130 m_colNormsUpdated(matrix.cols()),
131 m_colNormsDirect(matrix.cols()),
132 m_isInitialized(false),
134 {
135 compute(matrix.derived());
136 }
ColPivHouseholderQR & compute(const EigenBase< InputType > &matrix)

References Eigen::ColPivHouseholderQR< _MatrixType >::compute(), and Eigen::EigenBase< Derived >::derived().

+ Here is the call graph for this function:

◆ ColPivHouseholderQR() [4/4]

template<typename _MatrixType >
template<typename InputType >
Eigen::ColPivHouseholderQR< _MatrixType >::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&)
146 : m_qr(matrix.derived()),
147 m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
148 m_colsPermutation(PermIndexType(matrix.cols())),
149 m_colsTranspositions(matrix.cols()),
150 m_temp(matrix.cols()),
151 m_colNormsUpdated(matrix.cols()),
152 m_colNormsDirect(matrix.cols()),
153 m_isInitialized(false),
155 {
157 }
void computeInPlace()
Definition ColPivHouseholderQR.h:480

References Eigen::ColPivHouseholderQR< _MatrixType >::computeInPlace().

+ Here is the call graph for this function:

Member Function Documentation

◆ _solve_impl() [1/2]

template<typename _MatrixType >
template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void Eigen::ColPivHouseholderQR< _MatrixType >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const

◆ _solve_impl() [2/2]

template<typename _MatrixType >
template<typename RhsType , typename DstType >
void Eigen::ColPivHouseholderQR< _MatrixType >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const
586{
587 eigen_assert(rhs.rows() == rows());
588
589 const Index nonzero_pivots = nonzeroPivots();
590
591 if(nonzero_pivots == 0)
592 {
593 dst.setZero();
594 return;
595 }
596
597 typename RhsType::PlainObject c(rhs);
598
599 // Note that the matrix Q = H_0^* H_1^*... so its inverse is Q^* = (H_0 H_1 ...)^T
600 c.applyOnTheLeft(householderSequence(m_qr, m_hCoeffs)
601 .setLength(nonzero_pivots)
602 .transpose()
603 );
604
605 m_qr.topLeftCorner(nonzero_pivots, nonzero_pivots)
606 .template triangularView<Upper>()
607 .solveInPlace(c.topRows(nonzero_pivots));
608
609 for(Index i = 0; i < nonzero_pivots; ++i) dst.row(m_colsPermutation.indices().coeff(i)) = c.row(i);
610 for(Index i = nonzero_pivots; i < cols(); ++i) dst.row(m_colsPermutation.indices().coeff(i)).setZero();
611}
#define eigen_assert(x)
Definition Macros.h:579
Index nonzeroPivots() const
Definition ColPivHouseholderQR.h:394
const IndicesType & indices() const
Definition PermutationMatrix.h:388
HouseholderSequence< VectorsType, CoeffsType > householderSequence(const VectorsType &v, const CoeffsType &h)
Convenience function for constructing a Householder sequence.
Definition HouseholderSequence.h:451
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:33

References eigen_assert, and Eigen::householderSequence().

+ Here is the call graph for this function:

◆ absDeterminant()

template<typename MatrixType >
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType >::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
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 abs(m_qr.diagonal().prod());
454}
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half abs(const half &a)
Definition Half.h:445

References eigen_assert.

◆ check_template_parameters()

template<typename _MatrixType >
static void Eigen::ColPivHouseholderQR< _MatrixType >::check_template_parameters ( )
inlinestaticprotected
428 {
430 }
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition StaticAssert.h:184
MatrixType::Scalar Scalar
Definition ColPivHouseholderQR.h:59

References EIGEN_STATIC_ASSERT_NON_INTEGER.

◆ cols()

template<typename _MatrixType >
Index Eigen::ColPivHouseholderQR< _MatrixType >::cols ( ) const
inline
328{ return m_qr.cols(); }

References Eigen::ColPivHouseholderQR< _MatrixType >::m_qr.

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::cols(), Eigen::ColPivHouseholderQR< _MatrixType >::dimensionOfKernel(), Eigen::ColPivHouseholderQR< _MatrixType >::isInjective(), and Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixZ().

+ Here is the caller graph for this function:

◆ colsPermutation()

template<typename _MatrixType >
const PermutationType & Eigen::ColPivHouseholderQR< _MatrixType >::colsPermutation ( ) const
inline
Returns
a const reference to the column permutation matrix
215 {
216 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
217 return m_colsPermutation;
218 }

References eigen_assert, Eigen::ColPivHouseholderQR< _MatrixType >::m_colsPermutation, and Eigen::ColPivHouseholderQR< _MatrixType >::m_isInitialized.

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::colsPermutation().

+ Here is the caller graph for this function:

◆ compute() [1/2]

template<typename _MatrixType >
template<typename InputType >
ColPivHouseholderQR & Eigen::ColPivHouseholderQR< _MatrixType >::compute ( const EigenBase< InputType > &  matrix)

Referenced by Eigen::ColPivHouseholderQR< _MatrixType >::ColPivHouseholderQR(), and Eigen::CompleteOrthogonalDecomposition< _MatrixType >::compute().

+ Here is the caller graph for this function:

◆ compute() [2/2]

template<typename _MatrixType >
template<typename InputType >
ColPivHouseholderQR< MatrixType > & Eigen::ColPivHouseholderQR< _MatrixType >::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&)
473{
474 m_qr = matrix.derived();
476 return *this;
477}

References Eigen::EigenBase< Derived >::derived().

+ Here is the call graph for this function:

◆ computeInPlace()

template<typename MatrixType >
void Eigen::ColPivHouseholderQR< MatrixType >::computeInPlace
protected
481{
483
484 // the column permutation is stored as int indices, so just to be sure:
485 eigen_assert(m_qr.cols()<=NumTraits<int>::highest());
486
487 using std::abs;
488
489 Index rows = m_qr.rows();
490 Index cols = m_qr.cols();
491 Index size = m_qr.diagonalSize();
492
493 m_hCoeffs.resize(size);
494
495 m_temp.resize(cols);
496
497 m_colsTranspositions.resize(m_qr.cols());
498 Index number_of_transpositions = 0;
499
500 m_colNormsUpdated.resize(cols);
501 m_colNormsDirect.resize(cols);
502 for (Index k = 0; k < cols; ++k) {
503 // colNormsDirect(k) caches the most recent directly computed norm of
504 // column k.
505 m_colNormsDirect.coeffRef(k) = m_qr.col(k).norm();
506 m_colNormsUpdated.coeffRef(k) = m_colNormsDirect.coeffRef(k);
507 }
508
509 RealScalar threshold_helper = numext::abs2<RealScalar>(m_colNormsUpdated.maxCoeff() * NumTraits<RealScalar>::epsilon()) / RealScalar(rows);
510 RealScalar norm_downdate_threshold = numext::sqrt(NumTraits<RealScalar>::epsilon());
511
512 m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
514
515 for(Index k = 0; k < size; ++k)
516 {
517 // first, we look up in our table m_colNormsUpdated which column has the biggest norm
518 Index biggest_col_index;
519 RealScalar biggest_col_sq_norm = numext::abs2(m_colNormsUpdated.tail(cols-k).maxCoeff(&biggest_col_index));
520 biggest_col_index += k;
521
522 // Track the number of meaningful pivots but do not stop the decomposition to make
523 // sure that the initial matrix is properly reproduced. See bug 941.
524 if(m_nonzero_pivots==size && biggest_col_sq_norm < threshold_helper * RealScalar(rows-k))
526
527 // apply the transposition to the columns
528 m_colsTranspositions.coeffRef(k) = biggest_col_index;
529 if(k != biggest_col_index) {
530 m_qr.col(k).swap(m_qr.col(biggest_col_index));
531 std::swap(m_colNormsUpdated.coeffRef(k), m_colNormsUpdated.coeffRef(biggest_col_index));
532 std::swap(m_colNormsDirect.coeffRef(k), m_colNormsDirect.coeffRef(biggest_col_index));
533 ++number_of_transpositions;
534 }
535
536 // generate the householder vector, store it below the diagonal
537 RealScalar beta;
538 m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
539
540 // apply the householder transformation to the diagonal coefficient
541 m_qr.coeffRef(k,k) = beta;
542
543 // remember the maximum absolute value of diagonal coefficients
544 if(abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
545
546 // apply the householder transformation
547 m_qr.bottomRightCorner(rows-k, cols-k-1)
548 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
549
550 // update our table of norms of the columns
551 for (Index j = k + 1; j < cols; ++j) {
552 // The following implements the stable norm downgrade step discussed in
553 // http://www.netlib.org/lapack/lawnspdf/lawn176.pdf
554 // and used in LAPACK routines xGEQPF and xGEQP3.
555 // See lines 278-297 in http://www.netlib.org/lapack/explore-html/dc/df4/sgeqpf_8f_source.html
556 if (m_colNormsUpdated.coeffRef(j) != RealScalar(0)) {
557 RealScalar temp = abs(m_qr.coeffRef(k, j)) / m_colNormsUpdated.coeffRef(j);
558 temp = (RealScalar(1) + temp) * (RealScalar(1) - temp);
559 temp = temp < RealScalar(0) ? RealScalar(0) : temp;
560 RealScalar temp2 = temp * numext::abs2<RealScalar>(m_colNormsUpdated.coeffRef(j) /
561 m_colNormsDirect.coeffRef(j));
562 if (temp2 <= norm_downdate_threshold) {
563 // The updated norm has become too inaccurate so re-compute the column
564 // norm directly.
565 m_colNormsDirect.coeffRef(j) = m_qr.col(j).tail(rows - k - 1).norm();
566 m_colNormsUpdated.coeffRef(j) = m_colNormsDirect.coeffRef(j);
567 } else {
568 m_colNormsUpdated.coeffRef(j) *= numext::sqrt(temp);
569 }
570 }
571 }
572 }
573
575 for(PermIndexType k = 0; k < size/*m_nonzero_pivots*/; ++k)
577
578 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
579 m_isInitialized = true;
580}
Index m_det_pq
Definition ColPivHouseholderQR.h:444
MatrixType::RealScalar RealScalar
Definition ColPivHouseholderQR.h:60
static void check_template_parameters()
Definition ColPivHouseholderQR.h:427
RealScalar m_maxpivot
Definition ColPivHouseholderQR.h:442
Index m_nonzero_pivots
Definition ColPivHouseholderQR.h:443
Derived & applyTranspositionOnTheRight(Index i, Index j)
Definition PermutationMatrix.h:185
void setIdentity()
Definition PermutationMatrix.h:142
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE float sqrt(const float &x)
Definition MathFunctions.h:540
constexpr auto size(const C &c) -> decltype(c.size())
Definition span.hpp:183

References eigen_assert, and Eigen::numext::sqrt().

Referenced by Eigen::ColPivHouseholderQR< _MatrixType >::ColPivHouseholderQR().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ dimensionOfKernel()

template<typename _MatrixType >
Index Eigen::ColPivHouseholderQR< _MatrixType >::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&).
273 {
274 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
275 return cols() - rank();
276 }
Index rank() const
Definition ColPivHouseholderQR.h:255

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

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::dimensionOfKernel().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ hCoeffs()

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

For advanced uses only.

334{ return m_hCoeffs; }

References Eigen::ColPivHouseholderQR< _MatrixType >::m_hCoeffs.

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::hCoeffs().

+ Here is the caller graph for this function:

◆ 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())
635{
636 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
637 return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
638}
HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::type > HouseholderSequenceType
Definition ColPivHouseholderQR.h:68

References eigen_assert.

Referenced by Eigen::ColPivHouseholderQR< _MatrixType >::matrixQ(), and Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixQ().

+ Here is the caller graph for this function:

◆ info()

template<typename _MatrixType >
ComputationInfo Eigen::ColPivHouseholderQR< _MatrixType >::info ( ) const
inline

Reports whether the QR factorization was succesful.

Note
This function always returns Success. It is provided for compatibility with other factorization routines.
Returns
Success
412 {
413 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
414 return Success;
415 }
@ Success
Definition Constants.h:432

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

◆ inverse()

template<typename _MatrixType >
const Inverse< ColPivHouseholderQR > Eigen::ColPivHouseholderQR< _MatrixType >::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.
322 {
323 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
324 return Inverse<ColPivHouseholderQR>(*this);
325 }

References eigen_assert, and Eigen::ColPivHouseholderQR< _MatrixType >::m_isInitialized.

◆ isInjective()

template<typename _MatrixType >
bool Eigen::ColPivHouseholderQR< _MatrixType >::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&).
286 {
287 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
288 return rank() == cols();
289 }

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

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::isInjective(), and Eigen::ColPivHouseholderQR< _MatrixType >::isInvertible().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ isInvertible()

template<typename _MatrixType >
bool Eigen::ColPivHouseholderQR< _MatrixType >::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 isInjective() const
Definition ColPivHouseholderQR.h:285
bool isSurjective() const
Definition ColPivHouseholderQR.h:298

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

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::isInvertible().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ isSurjective()

template<typename _MatrixType >
bool Eigen::ColPivHouseholderQR< _MatrixType >::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&).
299 {
300 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
301 return rank() == rows();
302 }

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

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

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ logAbsDeterminant()

template<typename MatrixType >
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType >::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
absDeterminant(), MatrixBase::determinant()
458{
459 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
460 eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
461 return m_qr.diagonal().cwiseAbs().array().log().sum();
462}

References eigen_assert.

◆ matrixQ()

template<typename _MatrixType >
HouseholderSequenceType Eigen::ColPivHouseholderQR< _MatrixType >::matrixQ ( ) const
inline
183 {
184 return householderQ();
185 }
HouseholderSequenceType householderQ() const
Definition ColPivHouseholderQR.h:634

References Eigen::ColPivHouseholderQR< _MatrixType >::householderQ().

+ Here is the call graph for this function:

◆ matrixQR()

template<typename _MatrixType >
const MatrixType & Eigen::ColPivHouseholderQR< _MatrixType >::matrixQR ( ) const
inline
Returns
a reference to the matrix where the Householder QR decomposition is stored
190 {
191 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
192 return m_qr;
193 }

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

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixQTZ(), and Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixT().

+ Here is the caller graph for this function:

◆ matrixR()

template<typename _MatrixType >
const MatrixType & Eigen::ColPivHouseholderQR< _MatrixType >::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:204
For rank-deficient matrices, use
matrixR().topLeftCorner(rank(), rank()).template triangularView<Upper>()
205 {
206 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
207 return m_qr;
208 }

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

◆ maxPivot()

template<typename _MatrixType >
RealScalar Eigen::ColPivHouseholderQR< _MatrixType >::maxPivot ( ) const
inline
Returns
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of R.
403{ return m_maxpivot; }

References Eigen::ColPivHouseholderQR< _MatrixType >::m_maxpivot.

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::maxPivot().

+ Here is the caller graph for this function:

◆ nonzeroPivots()

template<typename _MatrixType >
Index Eigen::ColPivHouseholderQR< _MatrixType >::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()
395 {
396 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
397 return m_nonzero_pivots;
398 }

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

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::nonzeroPivots().

+ Here is the caller graph for this function:

◆ rank()

template<typename _MatrixType >
Index Eigen::ColPivHouseholderQR< _MatrixType >::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&).
256 {
257 using std::abs;
258 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
259 RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
260 Index result = 0;
261 for(Index i = 0; i < m_nonzero_pivots; ++i)
262 result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold);
263 return result;
264 }
RealScalar threshold() const
Definition ColPivHouseholderQR.h:378

References eigen_assert, Eigen::ColPivHouseholderQR< _MatrixType >::m_isInitialized, Eigen::ColPivHouseholderQR< _MatrixType >::m_maxpivot, Eigen::ColPivHouseholderQR< _MatrixType >::m_nonzero_pivots, Eigen::ColPivHouseholderQR< _MatrixType >::m_qr, and Eigen::ColPivHouseholderQR< _MatrixType >::threshold().

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

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ rows()

template<typename _MatrixType >
Index Eigen::ColPivHouseholderQR< _MatrixType >::rows ( ) const
inline
327{ return m_qr.rows(); }

References Eigen::ColPivHouseholderQR< _MatrixType >::m_qr.

Referenced by Eigen::ColPivHouseholderQR< _MatrixType >::isSurjective(), and Eigen::CompleteOrthogonalDecomposition< _MatrixType >::rows().

+ Here is the caller graph for this function:

◆ setThreshold() [1/2]

template<typename _MatrixType >
ColPivHouseholderQR & Eigen::ColPivHouseholderQR< _MatrixType >::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)

354 {
357 return *this;
358 }
RealScalar m_prescribedThreshold
Definition ColPivHouseholderQR.h:442

Referenced by Eigen::CompleteOrthogonalDecomposition< _MatrixType >::setThreshold(), and Eigen::CompleteOrthogonalDecomposition< _MatrixType >::setThreshold().

+ Here is the caller graph for this function:

◆ setThreshold() [2/2]

template<typename _MatrixType >
ColPivHouseholderQR & Eigen::ColPivHouseholderQR< _MatrixType >::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:352

See the documentation of setThreshold(const RealScalar&).

369 {
371 return *this;
372 }

References Eigen::ColPivHouseholderQR< _MatrixType >::m_usePrescribedThreshold.

◆ solve()

template<typename _MatrixType >
template<typename Rhs >
const Solve< ColPivHouseholderQR, Rhs > Eigen::ColPivHouseholderQR< _MatrixType >::solve ( const MatrixBase< Rhs > &  b) const
inline

This method finds a solution x to the equation Ax=b, where A is the matrix of which *this is the QR decomposition, if any exists.

Parameters
bthe right-hand-side of the equation to solve.
Returns
a solution.

\note_about_checking_solutions

\note_about_arbitrary_choice_of_solution

Example:

Output:

 
176 {
177 eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
178 return Solve<ColPivHouseholderQR, Rhs>(*this, b.derived());
179 }

References eigen_assert, and Eigen::ColPivHouseholderQR< _MatrixType >::m_isInitialized.

Referenced by igl::mvc().

+ Here is the caller graph for this function:

◆ threshold()

template<typename _MatrixType >
RealScalar Eigen::ColPivHouseholderQR< _MatrixType >::threshold ( ) const
inline

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

See the documentation of setThreshold(const RealScalar&).

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

References eigen_assert, Eigen::ColPivHouseholderQR< _MatrixType >::m_isInitialized, Eigen::ColPivHouseholderQR< _MatrixType >::m_prescribedThreshold, Eigen::ColPivHouseholderQR< _MatrixType >::m_qr, and Eigen::ColPivHouseholderQR< _MatrixType >::m_usePrescribedThreshold.

Referenced by Eigen::MatrixBase< Homogeneous< MatrixType, _Direction > >::colPivHouseholderQr(), Eigen::ColPivHouseholderQR< _MatrixType >::rank(), and Eigen::CompleteOrthogonalDecomposition< _MatrixType >::threshold().

+ Here is the caller graph for this function:

Friends And Related Symbol Documentation

◆ CompleteOrthogonalDecomposition< MatrixType >

template<typename _MatrixType >
friend class CompleteOrthogonalDecomposition< MatrixType >
friend

Member Data Documentation

◆ m_colNormsDirect

template<typename _MatrixType >
RealRowVectorType Eigen::ColPivHouseholderQR< _MatrixType >::m_colNormsDirect
protected

◆ m_colNormsUpdated

template<typename _MatrixType >
RealRowVectorType Eigen::ColPivHouseholderQR< _MatrixType >::m_colNormsUpdated
protected

◆ m_colsPermutation

template<typename _MatrixType >
PermutationType Eigen::ColPivHouseholderQR< _MatrixType >::m_colsPermutation
protected

◆ m_colsTranspositions

template<typename _MatrixType >
IntRowVectorType Eigen::ColPivHouseholderQR< _MatrixType >::m_colsTranspositions
protected

◆ m_det_pq

template<typename _MatrixType >
Index Eigen::ColPivHouseholderQR< _MatrixType >::m_det_pq
protected

◆ m_hCoeffs

template<typename _MatrixType >
HCoeffsType Eigen::ColPivHouseholderQR< _MatrixType >::m_hCoeffs
protected

◆ m_isInitialized

◆ m_maxpivot

template<typename _MatrixType >
RealScalar Eigen::ColPivHouseholderQR< _MatrixType >::m_maxpivot
protected

◆ m_nonzero_pivots

template<typename _MatrixType >
Index Eigen::ColPivHouseholderQR< _MatrixType >::m_nonzero_pivots
protected

◆ m_prescribedThreshold

◆ m_qr

◆ m_temp

template<typename _MatrixType >
RowVectorType Eigen::ColPivHouseholderQR< _MatrixType >::m_temp
protected

◆ m_usePrescribedThreshold


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