Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Private Types | Friends | List of all members
Eigen::CompleteOrthogonalDecomposition Class Reference

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

#include <ForwardDeclarations.h>

Public Types

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

Public Member Functions

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
 
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)
 
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)
 
RealScalar threshold () const
 
const HCoeffsTypezCoeffs () 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 ()
 

Static Protected Member Functions

static void check_template_parameters ()
 

Protected Attributes

ColPivHouseholderQR< MatrixTypem_cpqr
 
RowVectorType m_temp
 
HCoeffsType m_zCoeffs
 

Private Types

typedef PermutationType::Index PermIndexType
 

Friends

template<typename Derived >
struct internal::solve_assertion
 

Detailed Description

Complete orthogonal decomposition (COD) of a matrix.

Parameters
MatrixTypethe 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 276 of file ForwardDeclarations.h.

Member Typedef Documentation

◆ Base

Definition at line 55 of file CompleteOrthogonalDecomposition.h.

◆ HCoeffsType

Definition at line 65 of file CompleteOrthogonalDecomposition.h.

◆ HouseholderSequenceType

Definition at line 76 of file CompleteOrthogonalDecomposition.h.

◆ IntRowVectorType

Definition at line 69 of file CompleteOrthogonalDecomposition.h.

◆ MatrixType

Definition at line 54 of file CompleteOrthogonalDecomposition.h.

◆ PermIndexType

Definition at line 80 of file CompleteOrthogonalDecomposition.h.

◆ PermutationType

Definition at line 67 of file CompleteOrthogonalDecomposition.h.

◆ PlainObject

typedef MatrixType::PlainObject Eigen::CompleteOrthogonalDecomposition::PlainObject

Definition at line 77 of file CompleteOrthogonalDecomposition.h.

◆ RealRowVectorType

Definition at line 72 of file CompleteOrthogonalDecomposition.h.

◆ RowVectorType

Definition at line 70 of file CompleteOrthogonalDecomposition.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 61 of file CompleteOrthogonalDecomposition.h.

Constructor & Destructor Documentation

◆ CompleteOrthogonalDecomposition() [1/4]

Eigen::CompleteOrthogonalDecomposition::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 90 of file CompleteOrthogonalDecomposition.h.

◆ CompleteOrthogonalDecomposition() [2/4]

Eigen::CompleteOrthogonalDecomposition::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 98 of file CompleteOrthogonalDecomposition.h.

◆ CompleteOrthogonalDecomposition() [3/4]

template<typename InputType >
Eigen::CompleteOrthogonalDecomposition::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);
See also
compute()

Definition at line 118 of file CompleteOrthogonalDecomposition.h.

◆ CompleteOrthogonalDecomposition() [4/4]

template<typename InputType >
Eigen::CompleteOrthogonalDecomposition::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 133 of file CompleteOrthogonalDecomposition.h.

Member Function Documentation

◆ _check_solve_assertion()

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

Definition at line 385 of file CompleteOrthogonalDecomposition.h.

◆ _solve_impl()

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

Definition at line 534 of file CompleteOrthogonalDecomposition.h.

◆ _solve_impl_transposed()

template<bool Conjugate, typename RhsType , typename DstType >
void Eigen::CompleteOrthogonalDecomposition::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const

Definition at line 566 of file CompleteOrthogonalDecomposition.h.

◆ absDeterminant()

MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition::absDeterminant ( ) const
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
logAbsDeterminant(), MatrixBase::determinant()

Definition at line 412 of file CompleteOrthogonalDecomposition.h.

◆ applyZAdjointOnTheLeftInPlace()

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

Overwrites rhs with $ \mathbf{Z}^* * \mathbf{rhs} $.

Definition at line 511 of file CompleteOrthogonalDecomposition.h.

◆ applyZOnTheLeftInPlace()

template<bool Conjugate, typename Rhs >
void Eigen::CompleteOrthogonalDecomposition::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 489 of file CompleteOrthogonalDecomposition.h.

◆ check_template_parameters()

static void Eigen::CompleteOrthogonalDecomposition::check_template_parameters ( )
inlinestaticprotected

Definition at line 380 of file CompleteOrthogonalDecomposition.h.

◆ cols()

Index Eigen::CompleteOrthogonalDecomposition::cols ( ) const
inline

Definition at line 285 of file CompleteOrthogonalDecomposition.h.

◆ colsPermutation()

const PermutationType& Eigen::CompleteOrthogonalDecomposition::colsPermutation ( ) const
inline
Returns
a const reference to the column permutation matrix

Definition at line 194 of file CompleteOrthogonalDecomposition.h.

◆ compute()

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

Definition at line 186 of file CompleteOrthogonalDecomposition.h.

◆ computeInPlace()

void Eigen::CompleteOrthogonalDecomposition::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 430 of file CompleteOrthogonalDecomposition.h.

◆ dimensionOfKernel()

Index Eigen::CompleteOrthogonalDecomposition::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 244 of file CompleteOrthogonalDecomposition.h.

◆ hCoeffs()

const HCoeffsType& Eigen::CompleteOrthogonalDecomposition::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 292 of file CompleteOrthogonalDecomposition.h.

◆ householderQ()

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

Definition at line 619 of file CompleteOrthogonalDecomposition.h.

◆ info()

ComputationInfo Eigen::CompleteOrthogonalDecomposition::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 366 of file CompleteOrthogonalDecomposition.h.

◆ isInjective()

bool Eigen::CompleteOrthogonalDecomposition::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 253 of file CompleteOrthogonalDecomposition.h.

◆ isInvertible()

bool Eigen::CompleteOrthogonalDecomposition::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 271 of file CompleteOrthogonalDecomposition.h.

◆ isSurjective()

bool Eigen::CompleteOrthogonalDecomposition::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 262 of file CompleteOrthogonalDecomposition.h.

◆ logAbsDeterminant()

MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition::logAbsDeterminant ( ) const
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
absDeterminant(), MatrixBase::determinant()

Definition at line 418 of file CompleteOrthogonalDecomposition.h.

◆ matrixQ()

HouseholderSequenceType Eigen::CompleteOrthogonalDecomposition::matrixQ ( void  ) const
inline

Definition at line 157 of file CompleteOrthogonalDecomposition.h.

◆ matrixQTZ()

const MatrixType& Eigen::CompleteOrthogonalDecomposition::matrixQTZ ( ) const
inline
Returns
a reference to the matrix where the complete orthogonal decomposition is stored

Definition at line 170 of file CompleteOrthogonalDecomposition.h.

◆ matrixT()

const MatrixType& Eigen::CompleteOrthogonalDecomposition::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
matrixR().topLeftCorner(rank(), rank()).template triangularView<Upper>()

Definition at line 183 of file CompleteOrthogonalDecomposition.h.

◆ matrixZ()

MatrixType Eigen::CompleteOrthogonalDecomposition::matrixZ ( ) const
inline
Returns
the matrix Z.

Definition at line 161 of file CompleteOrthogonalDecomposition.h.

◆ maxPivot()

RealScalar Eigen::CompleteOrthogonalDecomposition::maxPivot ( ) const
inline
Returns
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of R.

Definition at line 356 of file CompleteOrthogonalDecomposition.h.

◆ nonzeroPivots()

Index Eigen::CompleteOrthogonalDecomposition::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 351 of file CompleteOrthogonalDecomposition.h.

◆ pseudoInverse()

const Inverse<CompleteOrthogonalDecomposition> Eigen::CompleteOrthogonalDecomposition::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 278 of file CompleteOrthogonalDecomposition.h.

◆ rank()

Index Eigen::CompleteOrthogonalDecomposition::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 235 of file CompleteOrthogonalDecomposition.h.

◆ rows()

Index Eigen::CompleteOrthogonalDecomposition::rows ( ) const
inline

Definition at line 284 of file CompleteOrthogonalDecomposition.h.

◆ setThreshold() [1/2]

CompleteOrthogonalDecomposition& Eigen::CompleteOrthogonalDecomposition::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 320 of file CompleteOrthogonalDecomposition.h.

◆ setThreshold() [2/2]

CompleteOrthogonalDecomposition& Eigen::CompleteOrthogonalDecomposition::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);

See the documentation of setThreshold(const RealScalar&).

Definition at line 333 of file CompleteOrthogonalDecomposition.h.

◆ threshold()

RealScalar Eigen::CompleteOrthogonalDecomposition::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 342 of file CompleteOrthogonalDecomposition.h.

◆ zCoeffs()

const HCoeffsType& Eigen::CompleteOrthogonalDecomposition::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 299 of file CompleteOrthogonalDecomposition.h.

Friends And Related Function Documentation

◆ internal::solve_assertion

template<typename Derived >
friend struct internal::solve_assertion
friend

Definition at line 58 of file CompleteOrthogonalDecomposition.h.

Member Data Documentation

◆ m_cpqr

ColPivHouseholderQR<MatrixType> Eigen::CompleteOrthogonalDecomposition::m_cpqr
protected

Definition at line 405 of file CompleteOrthogonalDecomposition.h.

◆ m_temp

RowVectorType Eigen::CompleteOrthogonalDecomposition::m_temp
protected

Definition at line 407 of file CompleteOrthogonalDecomposition.h.

◆ m_zCoeffs

HCoeffsType Eigen::CompleteOrthogonalDecomposition::m_zCoeffs
protected

Definition at line 406 of file CompleteOrthogonalDecomposition.h.


The documentation for this class was generated from the following files:
Eigen::CompleteOrthogonalDecomposition::cols
Index cols() const
Definition: CompleteOrthogonalDecomposition.h:285
Eigen::CompleteOrthogonalDecomposition::rank
Index rank() const
Definition: CompleteOrthogonalDecomposition.h:235
cod
void cod()
Definition: qr_colpivoting.cpp:17
Eigen::CompleteOrthogonalDecomposition::matrixT
const MatrixType & matrixT() const
Definition: CompleteOrthogonalDecomposition.h:183
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: gtsam/3rdparty/Eigen/blas/common.h:110
qr
HouseholderQR< MatrixXf > qr(A)
Eigen::Default
@ Default
Definition: Constants.h:362


gtsam
Author(s):
autogenerated on Thu Dec 19 2024 04:09:33