Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Private Types | List of all members
Eigen::CompleteOrthogonalDecomposition< _MatrixType > Class Template Reference

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

#include <ForwardDeclarations.h>

Public Types

enum  { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::type > HouseholderSequenceType
 
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 MatrixType::RealScalar RealScalar
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 
MatrixType::RealScalar absDeterminant () const
 
Index cols () const
 
const PermutationTypecolsPermutation () const
 
 CompleteOrthogonalDecomposition ()
 Default Constructor. More...
 
 CompleteOrthogonalDecomposition (Index rows, Index cols)
 Default Constructor with memory preallocation. 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...
 
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 succesful. 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
 

Protected Member Functions

template<typename Rhs >
void applyZAdjointOnTheLeftInPlace (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
 

Detailed Description

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

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 257 of file ForwardDeclarations.h.

Member Typedef Documentation

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

Definition at line 60 of file CompleteOrthogonalDecomposition.h.

Definition at line 71 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 64 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 50 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
typedef PermutationType::Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::PermIndexType
private

Definition at line 75 of file CompleteOrthogonalDecomposition.h.

Definition at line 62 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 72 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 67 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 58 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 65 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 57 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 59 of file CompleteOrthogonalDecomposition.h.

Member Enumeration Documentation

template<typename _MatrixType>
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 51 of file CompleteOrthogonalDecomposition.h.

Constructor & Destructor Documentation

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

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

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

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

Member Function Documentation

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

Definition at line 492 of file CompleteOrthogonalDecomposition.h.

template<typename MatrixType >
MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< MatrixType >::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 392 of file CompleteOrthogonalDecomposition.h.

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

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

Definition at line 469 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
static void Eigen::CompleteOrthogonalDecomposition< _MatrixType >::check_template_parameters ( )
inlinestaticprotected

Definition at line 374 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::cols ( void  ) const
inline

Definition at line 282 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 192 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 184 of file CompleteOrthogonalDecomposition.h.

template<typename MatrixType >
void Eigen::CompleteOrthogonalDecomposition< MatrixType >::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 410 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 242 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
const HCoeffsType& Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 289 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 546 of file CompleteOrthogonalDecomposition.h.

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

Reports whether the complete orthogonal decomposition was succesful.

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

Definition at line 363 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 251 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 269 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 260 of file CompleteOrthogonalDecomposition.h.

template<typename MatrixType >
MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< MatrixType >::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 398 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
HouseholderSequenceType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixQ ( void  ) const
inline

Definition at line 155 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 168 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
const MatrixType& Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 181 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 159 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 353 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 348 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
const Inverse<CompleteOrthogonalDecomposition> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 276 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 233 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::rows ( void  ) const
inline

Definition at line 281 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
CompleteOrthogonalDecomposition& Eigen::CompleteOrthogonalDecomposition< _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. 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 317 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
CompleteOrthogonalDecomposition& Eigen::CompleteOrthogonalDecomposition< _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);

See the documentation of setThreshold(const RealScalar&).

Definition at line 330 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
template<typename Rhs >
const Solve<CompleteOrthogonalDecomposition, Rhs> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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.

Definition at line 147 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
RealScalar Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 339 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
const HCoeffsType& Eigen::CompleteOrthogonalDecomposition< _MatrixType >::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 296 of file CompleteOrthogonalDecomposition.h.

Member Data Documentation

template<typename _MatrixType>
ColPivHouseholderQR<MatrixType> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::m_cpqr
protected

Definition at line 385 of file CompleteOrthogonalDecomposition.h.

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

Definition at line 387 of file CompleteOrthogonalDecomposition.h.

template<typename _MatrixType>
HCoeffsType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::m_zCoeffs
protected

Definition at line 386 of file CompleteOrthogonalDecomposition.h.


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


hebiros
Author(s): Xavier Artache , Matthew Tesch
autogenerated on Thu Sep 3 2020 04:09:59