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

Complete orthogonal decomposition (COD) of a matrix. More...

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

Public Types

enum  {
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef SolverBase< CompleteOrthogonalDecompositionBase
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef HouseholderSequence< MatrixType, internal::remove_all_t< typename HCoeffsType::ConjugateReturnType > > HouseholderSequenceType
 
typedef internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
 
typedef MatrixType_ MatrixType
 
typedef PermutationIndex_ PermutationIndex
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndexPermutationType
 
typedef MatrixType::PlainObject PlainObject
 
typedef internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
- Public Types inherited from Eigen::SolverBase< CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > >
enum  
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
typedef EigenBase< CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > > Base
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > > ConstTransposeReturnType
 
typedef internal::traits< CompleteOrthogonalDecomposition< 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
 
 CompleteOrthogonalDecomposition ()
 Default Constructor. More...
 
template<typename InputType >
 CompleteOrthogonalDecomposition (const EigenBase< InputType > &matrix)
 Constructs a complete orthogonal decomposition from a given matrix. More...
 
template<typename InputType >
 CompleteOrthogonalDecomposition (EigenBase< InputType > &matrix)
 Constructs a complete orthogonal decomposition from a given matrix. More...
 
 CompleteOrthogonalDecomposition (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
template<typename InputType >
CompleteOrthogonalDecompositioncompute (const EigenBase< InputType > &matrix)
 
MatrixType::Scalar determinant () const
 
Index dimensionOfKernel () const
 
const HCoeffsTypehCoeffs () const
 
HouseholderSequenceType householderQ (void) const
 
ComputationInfo info () const
 Reports whether the complete orthogonal decomposition was successful. More...
 
bool isInjective () const
 
bool isInvertible () const
 
bool isSurjective () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
HouseholderSequenceType matrixQ (void) const
 
const MatrixTypematrixQTZ () const
 
const MatrixTypematrixT () const
 
MatrixType matrixZ () const
 
RealScalar maxPivot () const
 
Index nonzeroPivots () const
 
const Inverse< CompleteOrthogonalDecompositionpseudoInverse () const
 
Index rank () const
 
Index rows () const
 
CompleteOrthogonalDecompositionsetThreshold (const RealScalar &threshold)
 
CompleteOrthogonalDecompositionsetThreshold (Default_t)
 
template<typename Rhs >
const Solve< CompleteOrthogonalDecomposition, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
const HCoeffsTypezCoeffs () const
 
- Public Member Functions inherited from Eigen::SolverBase< CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > >
const AdjointReturnType adjoint () const
 
CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > & derived ()
 
const CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > & derived () const
 
const Solve< CompleteOrthogonalDecomposition< 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

template<bool Transpose_, typename Rhs >
void _check_solve_assertion (const Rhs &b) const
 
template<typename Rhs >
void applyZAdjointOnTheLeftInPlace (Rhs &rhs) const
 
template<bool Conjugate, typename Rhs >
void applyZOnTheLeftInPlace (Rhs &rhs) const
 
void computeInPlace ()
 
- Protected Member Functions inherited from Eigen::SolverBase< CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > >
void _check_solve_assertion (const Rhs &b) const
 

Protected Attributes

ColPivHouseholderQR< MatrixType, PermutationIndexm_cpqr
 
RowVectorType m_temp
 
HCoeffsType m_zCoeffs
 

Detailed Description

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

Complete orthogonal decomposition (COD) of a matrix.

Template Parameters
MatrixType_the type of the matrix of which we are computing the COD.

This class performs a rank-revealing complete orthogonal decomposition of a matrix A into matrices P, Q, T, and Z such that

\[ \mathbf{A} \, \mathbf{P} = \mathbf{Q} \, \begin{bmatrix} \mathbf{T} & \mathbf{0} \\ \mathbf{0} & \mathbf{0} \end{bmatrix} \, \mathbf{Z} \]

by using Householder transformations. Here, P is a permutation matrix, Q and Z are unitary matrices and T an upper triangular matrix of size rank-by-rank. A may be rank deficient.

This class supports the inplace decomposition mechanism.

See also
MatrixBase::completeOrthogonalDecomposition()

Definition at line 52 of file CompleteOrthogonalDecomposition.h.

Member Typedef Documentation

◆ Base

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

Definition at line 57 of file CompleteOrthogonalDecomposition.h.

◆ HCoeffsType

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

Definition at line 67 of file CompleteOrthogonalDecomposition.h.

◆ HouseholderSequenceType

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

Definition at line 78 of file CompleteOrthogonalDecomposition.h.

◆ IntRowVectorType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::plain_row_type<MatrixType, Index>::type Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::IntRowVectorType

Definition at line 71 of file CompleteOrthogonalDecomposition.h.

◆ MatrixType

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

Definition at line 56 of file CompleteOrthogonalDecomposition.h.

◆ PermutationIndex

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

Definition at line 61 of file CompleteOrthogonalDecomposition.h.

◆ PermutationType

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

Definition at line 69 of file CompleteOrthogonalDecomposition.h.

◆ PlainObject

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

Definition at line 79 of file CompleteOrthogonalDecomposition.h.

◆ RealRowVectorType

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

Definition at line 74 of file CompleteOrthogonalDecomposition.h.

◆ RowVectorType

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

Definition at line 72 of file CompleteOrthogonalDecomposition.h.

Member Enumeration Documentation

◆ anonymous enum

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

Definition at line 63 of file CompleteOrthogonalDecomposition.h.

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

Constructor & Destructor Documentation

◆ CompleteOrthogonalDecomposition() [1/4]

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

Default Constructor.

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

Definition at line 89 of file CompleteOrthogonalDecomposition.h.

◆ CompleteOrthogonalDecomposition() [2/4]

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

Definition at line 97 of file CompleteOrthogonalDecomposition.h.

◆ CompleteOrthogonalDecomposition() [3/4]

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

Constructs a complete orthogonal decomposition from a given matrix.

This constructor computes the complete orthogonal decomposition of the matrix matrix by calling the method compute(). The default threshold for rank determination will be used. It is a short cut for:

CompleteOrthogonalDecomposition<MatrixType> cod(matrix.rows(),
matrix.cols());
cod.setThreshold(Default);
cod.compute(matrix);
@ Default
Definition: Constants.h:364
See also
compute()

Definition at line 117 of file CompleteOrthogonalDecomposition.h.

118  : m_cpqr(matrix.rows(), matrix.cols()),
119  m_zCoeffs((std::min)(matrix.rows(), matrix.cols())),
120  m_temp(matrix.cols())
121  {
122  compute(matrix.derived());
123  }
CompleteOrthogonalDecomposition & compute(const EigenBase< InputType > &matrix)

◆ CompleteOrthogonalDecomposition() [4/4]

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename InputType >
Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::CompleteOrthogonalDecomposition ( EigenBase< InputType > &  matrix)
inlineexplicit

Constructs a complete orthogonal decomposition from a given matrix.

This overloaded constructor is provided for inplace decomposition when MatrixType is a Eigen::Ref.

See also
CompleteOrthogonalDecomposition(const EigenBase&)

Definition at line 132 of file CompleteOrthogonalDecomposition.h.

133  : m_cpqr(matrix.derived()),
134  m_zCoeffs((std::min)(matrix.rows(), matrix.cols())),
135  m_temp(matrix.cols())
136  {
137  computeInPlace();
138  }

Member Function Documentation

◆ _check_solve_assertion()

template<typename MatrixType_ , typename PermutationIndex_ >
template<bool Transpose_, typename Rhs >
void Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::_check_solve_assertion ( const Rhs &  b) const
inlineprotected

Definition at line 397 of file CompleteOrthogonalDecomposition.h.

397  {
399  eigen_assert(m_cpqr.m_isInitialized && "CompleteOrthogonalDecomposition is not initialized.");
400  eigen_assert((Transpose_?derived().cols():derived().rows())==b.rows() && "CompleteOrthogonalDecomposition::solve(): invalid number of rows of the right hand side matrix b");
401  }
Array< int, 3, 1 > b
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
Definition: Macros.h:914
#define eigen_assert(x)
Definition: Macros.h:902
CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > & derived()
Definition: EigenBase.h:48

◆ absDeterminant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >::absDeterminant
Returns
the absolute value of the determinant of the matrix of which *this is the complete orthogonal decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the complete orthogonal 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 430 of file CompleteOrthogonalDecomposition.h.

430  {
431  return m_cpqr.absDeterminant();
432 }
MatrixType::RealScalar absDeterminant() const

◆ applyZAdjointOnTheLeftInPlace()

template<typename MatrixType , typename PermutationIndex >
template<typename Rhs >
void Eigen::CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >::applyZAdjointOnTheLeftInPlace ( Rhs &  rhs) const
protected

Overwrites rhs with \( \mathbf{Z}^* * \mathbf{rhs} \).

Definition at line 526 of file CompleteOrthogonalDecomposition.h.

527  {
528  const Index cols = this->cols();
529  const Index nrhs = rhs.cols();
530  const Index rank = this->rank();
531  Matrix<typename Rhs::Scalar, Dynamic, 1> temp((std::max)(cols, nrhs));
532  for (Index k = 0; k < rank; ++k) {
533  if (k != rank - 1) {
534  rhs.row(k).swap(rhs.row(rank - 1));
535  }
536  rhs.middleRows(rank - 1, cols - rank + 1)
537  .applyHouseholderOnTheLeft(
538  matrixQTZ().row(k).tail(cols - rank).adjoint(), zCoeffs()(k),
539  &temp(0));
540  if (k != rank - 1) {
541  rhs.row(k).swap(rhs.row(rank - 1));
542  }
543  }
544 }
RowXpr row(Index i)
This is the const version of row(). *‍/.
FixedSegmentReturnType<... >::Type tail(NType n)
bfloat16() max(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:690
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41

◆ applyZOnTheLeftInPlace()

template<typename MatrixType , typename PermutationIndex >
template<bool Conjugate, typename Rhs >
void Eigen::CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >::applyZOnTheLeftInPlace ( Rhs &  rhs) const
protected

Overwrites rhs with \( \mathbf{Z} * \mathbf{rhs} \) or \( \mathbf{\overline Z} * \mathbf{rhs} \) if Conjugate is set to true.

Definition at line 504 of file CompleteOrthogonalDecomposition.h.

505  {
506  const Index cols = this->cols();
507  const Index nrhs = rhs.cols();
508  const Index rank = this->rank();
509  Matrix<typename Rhs::Scalar, Dynamic, 1> temp((std::max)(cols, nrhs));
510  for (Index k = rank-1; k >= 0; --k) {
511  if (k != rank - 1) {
512  rhs.row(k).swap(rhs.row(rank - 1));
513  }
514  rhs.middleRows(rank - 1, cols - rank + 1)
515  .applyHouseholderOnTheLeft(
516  matrixQTZ().row(k).tail(cols - rank).transpose().template conjugateIf<!Conjugate>(), zCoeffs().template conjugateIf<Conjugate>()(k),
517  &temp(0));
518  if (k != rank - 1) {
519  rhs.row(k).swap(rhs.row(rank - 1));
520  }
521  }
522 }

◆ cols()

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

Definition at line 299 of file CompleteOrthogonalDecomposition.h.

299 { return m_cpqr.cols(); }

◆ colsPermutation()

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

Definition at line 193 of file CompleteOrthogonalDecomposition.h.

193  {
194  return m_cpqr.colsPermutation();
195  }
const PermutationType & colsPermutation() const

◆ compute()

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

Definition at line 185 of file CompleteOrthogonalDecomposition.h.

185  {
186  // Compute the column pivoted QR factorization A P = Q R.
187  m_cpqr.compute(matrix);
188  computeInPlace();
189  return *this;
190  }
ColPivHouseholderQR & compute(const EigenBase< InputType > &matrix)

◆ computeInPlace()

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

Performs the complete orthogonal decomposition of the given matrix matrix. The result of the factorization is stored into *this, and a reference to *this is returned.

See also
class CompleteOrthogonalDecomposition, CompleteOrthogonalDecomposition(const MatrixType&)

Definition at line 448 of file CompleteOrthogonalDecomposition.h.

449 {
450  eigen_assert(m_cpqr.cols() <= NumTraits<PermutationIndex>::highest());
451 
452  const Index rank = m_cpqr.rank();
453  const Index cols = m_cpqr.cols();
454  const Index rows = m_cpqr.rows();
455  m_zCoeffs.resize((std::min)(rows, cols));
456  m_temp.resize(cols);
457 
458  if (rank < cols) {
459  // We have reduced the (permuted) matrix to the form
460  // [R11 R12]
461  // [ 0 R22]
462  // where R11 is r-by-r (r = rank) upper triangular, R12 is
463  // r-by-(n-r), and R22 is empty or the norm of R22 is negligible.
464  // We now compute the complete orthogonal decomposition by applying
465  // Householder transformations from the right to the upper trapezoidal
466  // matrix X = [R11 R12] to zero out R12 and obtain the factorization
467  // [R11 R12] = [T11 0] * Z, where T11 is r-by-r upper triangular and
468  // Z = Z(0) * Z(1) ... Z(r-1) is an n-by-n orthogonal matrix.
469  // We store the data representing Z in R12 and m_zCoeffs.
470  for (Index k = rank - 1; k >= 0; --k) {
471  if (k != rank - 1) {
472  // Given the API for Householder reflectors, it is more convenient if
473  // we swap the leading parts of columns k and r-1 (zero-based) to form
474  // the matrix X_k = [X(0:k, k), X(0:k, r:n)]
475  m_cpqr.m_qr.col(k).head(k + 1).swap(
476  m_cpqr.m_qr.col(rank - 1).head(k + 1));
477  }
478  // Construct Householder reflector Z(k) to zero out the last row of X_k,
479  // i.e. choose Z(k) such that
480  // [X(k, k), X(k, r:n)] * Z(k) = [beta, 0, .., 0].
481  RealScalar beta;
482  m_cpqr.m_qr.row(k)
483  .tail(cols - rank + 1)
484  .makeHouseholderInPlace(m_zCoeffs(k), beta);
485  m_cpqr.m_qr(k, rank - 1) = beta;
486  if (k > 0) {
487  // Apply Z(k) to the first k rows of X_k
488  m_cpqr.m_qr.topRightCorner(k, cols - rank + 1)
489  .applyHouseholderOnTheRight(
490  m_cpqr.m_qr.row(k).tail(cols - rank).adjoint(), m_zCoeffs(k),
491  &m_temp(0));
492  }
493  if (k != rank - 1) {
494  // Swap X(0:k,k) back to its proper location.
495  m_cpqr.m_qr.col(k).head(k + 1).swap(
496  m_cpqr.m_qr.col(rank - 1).head(k + 1));
497  }
498  }
499  }
500 }

◆ determinant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::Scalar Eigen::CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >::determinant
Returns
the determinant of the matrix of which *this is the complete orthogonal decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the complete orthogonal 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 424 of file CompleteOrthogonalDecomposition.h.

424  {
425  return m_cpqr.determinant();
426 }
MatrixType::Scalar determinant() const

◆ dimensionOfKernel()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::dimensionOfKernel ( ) const
inline
Returns
the dimension of the kernel of the matrix of which *this is the complete orthogonal 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 258 of file CompleteOrthogonalDecomposition.h.

258 { return m_cpqr.dimensionOfKernel(); }

◆ hCoeffs()

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

306 { return m_cpqr.hCoeffs(); }
const HCoeffsType & hCoeffs() const

◆ householderQ()

template<typename MatrixType , typename PermutationIndex >
CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >::HouseholderSequenceType Eigen::CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >::householderQ ( void  ) const
Returns
the matrix Q as a sequence of householder transformations

Definition at line 634 of file CompleteOrthogonalDecomposition.h.

634  {
635  return m_cpqr.householderQ();
636 }
HouseholderSequenceType householderQ() const

◆ info()

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

Reports whether the complete orthogonal decomposition was successful.

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

Definition at line 380 of file CompleteOrthogonalDecomposition.h.

380  {
381  eigen_assert(m_cpqr.m_isInitialized && "Decomposition is not initialized.");
382  return Success;
383  }
@ Success
Definition: Constants.h:446

◆ isInjective()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::isInjective ( ) const
inline
Returns
true if the matrix of which *this is the 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 267 of file CompleteOrthogonalDecomposition.h.

267 { return m_cpqr.isInjective(); }

◆ isInvertible()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::isInvertible ( ) const
inline
Returns
true if the matrix of which *this is the complete orthogonal 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 285 of file CompleteOrthogonalDecomposition.h.

285 { return m_cpqr.isInvertible(); }

◆ isSurjective()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::isSurjective ( ) const
inline
Returns
true if the matrix of which *this is the 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 276 of file CompleteOrthogonalDecomposition.h.

276 { return m_cpqr.isSurjective(); }

◆ logAbsDeterminant()

template<typename MatrixType , typename PermutationIndex >
MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< MatrixType, PermutationIndex >::logAbsDeterminant
Returns
the natural log of the absolute value of the determinant of the matrix of which *this is the complete orthogonal decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the complete orthogonal 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 436 of file CompleteOrthogonalDecomposition.h.

436  {
437  return m_cpqr.logAbsDeterminant();
438 }
MatrixType::RealScalar logAbsDeterminant() const

◆ matrixQ()

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

Definition at line 156 of file CompleteOrthogonalDecomposition.h.

156 { return m_cpqr.householderQ(); }

◆ matrixQTZ()

template<typename MatrixType_ , typename PermutationIndex_ >
const MatrixType& Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::matrixQTZ ( ) const
inline
Returns
a reference to the matrix where the complete orthogonal decomposition is stored

Definition at line 169 of file CompleteOrthogonalDecomposition.h.

169 { return m_cpqr.matrixQR(); }
const MatrixType & matrixQR() const

◆ matrixT()

template<typename MatrixType_ , typename PermutationIndex_ >
const MatrixType& Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::matrixT ( ) const
inline
Returns
a reference to the matrix where the complete orthogonal decomposition is stored.
Warning
The strict lower part and
cols() - rank()
right columns of this matrix contains internal values. Only the upper triangular part should be referenced. To get it, use
matrixT().template triangularView<Upper>()
For rank-deficient matrices, use
matrixT().topLeftCorner(rank(), rank()).template triangularView<Upper>()

Definition at line 182 of file CompleteOrthogonalDecomposition.h.

182 { return m_cpqr.matrixQR(); }

◆ matrixZ()

template<typename MatrixType_ , typename PermutationIndex_ >
MatrixType Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::matrixZ ( ) const
inline
Returns
the matrix Z.

Definition at line 160 of file CompleteOrthogonalDecomposition.h.

160  {
161  MatrixType Z = MatrixType::Identity(m_cpqr.cols(), m_cpqr.cols());
162  applyZOnTheLeftInPlace<false>(Z);
163  return Z;
164  }
Matrix< float, 1, Dynamic > MatrixType

◆ maxPivot()

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

Definition at line 370 of file CompleteOrthogonalDecomposition.h.

370 { return m_cpqr.maxPivot(); }

◆ nonzeroPivots()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::nonzeroPivots ( ) const
inline
Returns
the number of nonzero pivots in the complete orthogonal 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 365 of file CompleteOrthogonalDecomposition.h.

365 { return m_cpqr.nonzeroPivots(); }

◆ pseudoInverse()

template<typename MatrixType_ , typename PermutationIndex_ >
const Inverse<CompleteOrthogonalDecomposition> Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::pseudoInverse ( ) const
inline
Returns
the pseudo-inverse of the matrix of which *this is the complete orthogonal decomposition.
Warning
: Do not compute this->pseudoInverse()*rhs to solve a linear systems. It is more efficient and numerically stable to call this->solve(rhs).

Definition at line 292 of file CompleteOrthogonalDecomposition.h.

293  {
294  eigen_assert(m_cpqr.m_isInitialized && "CompleteOrthogonalDecomposition is not initialized.");
295  return Inverse<CompleteOrthogonalDecomposition>(*this);
296  }

◆ rank()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::rank ( ) const
inline
Returns
the rank of the matrix of which *this is the complete orthogonal 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 249 of file CompleteOrthogonalDecomposition.h.

249 { return m_cpqr.rank(); }

◆ rows()

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

Definition at line 298 of file CompleteOrthogonalDecomposition.h.

298 { return m_cpqr.rows(); }

◆ setThreshold() [1/2]

template<typename MatrixType_ , typename PermutationIndex_ >
CompleteOrthogonalDecomposition& Eigen::CompleteOrthogonalDecomposition< 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. Most be called before calling compute().

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 334 of file CompleteOrthogonalDecomposition.h.

334  {
336  return *this;
337  }
ColPivHouseholderQR & setThreshold(const RealScalar &threshold)

◆ setThreshold() [2/2]

template<typename MatrixType_ , typename PermutationIndex_ >
CompleteOrthogonalDecomposition& Eigen::CompleteOrthogonalDecomposition< 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);
HouseholderQR< MatrixXf > qr(A)

See the documentation of setThreshold(const RealScalar&).

Definition at line 347 of file CompleteOrthogonalDecomposition.h.

347  {
349  return *this;
350  }

◆ solve()

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename Rhs >
const Solve<CompleteOrthogonalDecomposition, Rhs> Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::solve ( const MatrixBase< Rhs > &  b) const
inline

This method computes the minimum-norm solution X to a least squares problem

\[\mathrm{minimize} \|A X - B\|, \]

where A is the matrix of which *this is the complete orthogonal decomposition.

Parameters
bthe right-hand sides of the problem to solve.
Returns
a solution.

◆ threshold()

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::CompleteOrthogonalDecomposition< 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 356 of file CompleteOrthogonalDecomposition.h.

356 { return m_cpqr.threshold(); }

◆ zCoeffs()

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

For advanced uses only.

Definition at line 313 of file CompleteOrthogonalDecomposition.h.

313 { return m_zCoeffs; }

Member Data Documentation

◆ m_cpqr

template<typename MatrixType_ , typename PermutationIndex_ >
ColPivHouseholderQR<MatrixType, PermutationIndex> Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::m_cpqr
protected

Definition at line 417 of file CompleteOrthogonalDecomposition.h.

◆ m_temp

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

Definition at line 419 of file CompleteOrthogonalDecomposition.h.

◆ m_zCoeffs

template<typename MatrixType_ , typename PermutationIndex_ >
HCoeffsType Eigen::CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >::m_zCoeffs
protected

Definition at line 418 of file CompleteOrthogonalDecomposition.h.


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