MOAB
4.9.3pre
|
Complete orthogonal decomposition (COD) of a matrix. More...
#include <CompleteOrthogonalDecomposition.h>
Public Types | |
enum | { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime } |
typedef _MatrixType | MatrixType |
typedef MatrixType::Scalar | Scalar |
typedef MatrixType::RealScalar | RealScalar |
typedef MatrixType::StorageIndex | StorageIndex |
typedef Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime > | MatrixQType |
typedef internal::plain_diag_type < MatrixType >::type | HCoeffsType |
typedef PermutationMatrix < ColsAtCompileTime, MaxColsAtCompileTime > | PermutationType |
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 | |
CompleteOrthogonalDecomposition () | |
Default Constructor. | |
CompleteOrthogonalDecomposition (Index rows, Index cols) | |
Default Constructor with memory preallocation. | |
template<typename InputType > | |
CompleteOrthogonalDecomposition (const EigenBase< InputType > &matrix) | |
Constructs a complete orthogonal decomposition from a given matrix. | |
template<typename Rhs > | |
const Solve < CompleteOrthogonalDecomposition, Rhs > | solve (const MatrixBase< Rhs > &b) const |
HouseholderSequenceType | householderQ (void) const |
HouseholderSequenceType | matrixQ (void) const |
MatrixType | matrixZ () const |
const MatrixType & | matrixQTZ () const |
const MatrixType & | matrixT () const |
template<typename InputType > | |
CompleteOrthogonalDecomposition & | compute (const EigenBase< InputType > &matrix) |
const PermutationType & | colsPermutation () 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 < CompleteOrthogonalDecomposition > | pseudoInverse () const |
Index | rows () const |
Index | cols () const |
const HCoeffsType & | hCoeffs () const |
const HCoeffsType & | zCoeffs () const |
CompleteOrthogonalDecomposition & | setThreshold (const RealScalar &threshold) |
CompleteOrthogonalDecomposition & | setThreshold (Default_t) |
RealScalar | threshold () const |
Index | nonzeroPivots () const |
RealScalar | maxPivot () const |
ComputationInfo | info () const |
Reports whether the complete orthogonal decomposition was succesful. | |
template<typename RhsType , typename DstType > | |
EIGEN_DEVICE_FUNC void | _solve_impl (const RhsType &rhs, DstType &dst) const |
Protected Member Functions | |
template<typename Rhs > | |
void | applyZAdjointOnTheLeftInPlace (Rhs &rhs) const |
Static Protected Member Functions | |
static void | check_template_parameters () |
Protected Attributes | |
ColPivHouseholderQR< MatrixType > | m_cpqr |
HCoeffsType | m_zCoeffs |
RowVectorType | m_temp |
Private Types | |
typedef PermutationType::Index | PermIndexType |
Complete orthogonal decomposition (COD) of a matrix.
MatrixType | the type of the matrix of which we are computing the COD. |
This class performs a rank-revealing complete ortogonal decomposition of a matrix A into matrices P, Q, T, and Z such that
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.
Definition at line 45 of file CompleteOrthogonalDecomposition.h.
typedef internal::plain_diag_type<MatrixType>::type Eigen::CompleteOrthogonalDecomposition< _MatrixType >::HCoeffsType |
Definition at line 61 of file CompleteOrthogonalDecomposition.h.
typedef HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType>::type> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::HouseholderSequenceType |
Definition at line 72 of file CompleteOrthogonalDecomposition.h.
typedef internal::plain_row_type<MatrixType, Index>::type Eigen::CompleteOrthogonalDecomposition< _MatrixType >::IntRowVectorType |
Definition at line 65 of file CompleteOrthogonalDecomposition.h.
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::MatrixQType |
Definition at line 60 of file CompleteOrthogonalDecomposition.h.
typedef _MatrixType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::MatrixType |
Definition at line 47 of file CompleteOrthogonalDecomposition.h.
typedef PermutationType::Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::PermIndexType [private] |
Definition at line 76 of file CompleteOrthogonalDecomposition.h.
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::PermutationType |
Definition at line 63 of file CompleteOrthogonalDecomposition.h.
typedef MatrixType::PlainObject Eigen::CompleteOrthogonalDecomposition< _MatrixType >::PlainObject |
Definition at line 73 of file CompleteOrthogonalDecomposition.h.
typedef internal::plain_row_type<MatrixType, RealScalar>::type Eigen::CompleteOrthogonalDecomposition< _MatrixType >::RealRowVectorType |
Definition at line 68 of file CompleteOrthogonalDecomposition.h.
typedef MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< _MatrixType >::RealScalar |
Definition at line 56 of file CompleteOrthogonalDecomposition.h.
typedef internal::plain_row_type<MatrixType>::type Eigen::CompleteOrthogonalDecomposition< _MatrixType >::RowVectorType |
Definition at line 66 of file CompleteOrthogonalDecomposition.h.
typedef MatrixType::Scalar Eigen::CompleteOrthogonalDecomposition< _MatrixType >::Scalar |
Definition at line 55 of file CompleteOrthogonalDecomposition.h.
typedef MatrixType::StorageIndex Eigen::CompleteOrthogonalDecomposition< _MatrixType >::StorageIndex |
Definition at line 57 of file CompleteOrthogonalDecomposition.h.
anonymous enum |
Definition at line 48 of file CompleteOrthogonalDecomposition.h.
{ RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
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 86 of file CompleteOrthogonalDecomposition.h.
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.
Definition at line 94 of file CompleteOrthogonalDecomposition.h.
Eigen::CompleteOrthogonalDecomposition< _MatrixType >::CompleteOrthogonalDecomposition | ( | const EigenBase< InputType > & | matrix | ) | [inline, explicit] |
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);
Definition at line 114 of file CompleteOrthogonalDecomposition.h.
void Eigen::CompleteOrthogonalDecomposition< _MatrixType >::_solve_impl | ( | const RhsType & | rhs, |
DstType & | dst | ||
) | const |
Definition at line 470 of file CompleteOrthogonalDecomposition.h.
{ eigen_assert(rhs.rows() == this->rows()); const Index rank = this->rank(); if (rank == 0) { dst.setZero(); return; } // Compute c = Q^* * rhs // Note that the matrix Q = H_0^* H_1^*... so its inverse is // Q^* = (H_0 H_1 ...)^T typename RhsType::PlainObject c(rhs); c.applyOnTheLeft( householderSequence(matrixQTZ(), hCoeffs()).setLength(rank).transpose()); // Solve T z = c(1:rank, :) dst.topRows(rank) = matrixT() .topLeftCorner(rank, rank) .template triangularView<Upper>() .solve(c.topRows(rank)); const Index cols = this->cols(); if (rank < cols) { // Compute y = Z^* * [ z ] // [ 0 ] dst.bottomRows(cols - rank).setZero(); applyZAdjointOnTheLeftInPlace(dst); } // Undo permutation to get x = P^{-1} * y. dst = colsPermutation() * dst; }
MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< MatrixType >::absDeterminant | ( | ) | const |
Definition at line 369 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.absDeterminant(); }
void Eigen::CompleteOrthogonalDecomposition< MatrixType >::applyZAdjointOnTheLeftInPlace | ( | Rhs & | rhs | ) | const [protected] |
Overwrites rhs with .
Definition at line 447 of file CompleteOrthogonalDecomposition.h.
{ const Index cols = this->cols(); const Index nrhs = rhs.cols(); const Index rank = this->rank(); Matrix<typename MatrixType::Scalar, Dynamic, 1> temp((std::max)(cols, nrhs)); for (Index k = 0; k < rank; ++k) { if (k != rank - 1) { rhs.row(k).swap(rhs.row(rank - 1)); } rhs.middleRows(rank - 1, cols - rank + 1) .applyHouseholderOnTheLeft( matrixQTZ().row(k).tail(cols - rank).adjoint(), zCoeffs()(k), &temp(0)); if (k != rank - 1) { rhs.row(k).swap(rhs.row(rank - 1)); } } }
static void Eigen::CompleteOrthogonalDecomposition< _MatrixType >::check_template_parameters | ( | ) | [inline, static, protected] |
Definition at line 353 of file CompleteOrthogonalDecomposition.h.
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::cols | ( | void | ) | const [inline] |
Definition at line 261 of file CompleteOrthogonalDecomposition.h.
const PermutationType& Eigen::CompleteOrthogonalDecomposition< _MatrixType >::colsPermutation | ( | ) | const [inline] |
Definition at line 171 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.colsPermutation(); }
CompleteOrthogonalDecomposition< MatrixType > & Eigen::CompleteOrthogonalDecomposition< MatrixType >::compute | ( | const EigenBase< InputType > & | matrix | ) |
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.
Definition at line 389 of file CompleteOrthogonalDecomposition.h.
{ check_template_parameters(); // the column permutation is stored as int indices, so just to be sure: eigen_assert(matrix.cols() <= NumTraits<int>::highest()); // Compute the column pivoted QR factorization A P = Q R. m_cpqr.compute(matrix); const Index rank = m_cpqr.rank(); const Index cols = matrix.cols(); if (rank < cols) { // We have reduced the (permuted) matrix to the form // [R11 R12] // [ 0 R22] // where R11 is r-by-r (r = rank) upper triangular, R12 is // r-by-(n-r), and R22 is empty or the norm of R22 is negligible. // We now compute the complete orthogonal decomposition by applying // Householder transformations from the right to the upper trapezoidal // matrix X = [R11 R12] to zero out R12 and obtain the factorization // [R11 R12] = [T11 0] * Z, where T11 is r-by-r upper triangular and // Z = Z(0) * Z(1) ... Z(r-1) is an n-by-n orthogonal matrix. // We store the data representing Z in R12 and m_zCoeffs. for (Index k = rank - 1; k >= 0; --k) { if (k != rank - 1) { // Given the API for Householder reflectors, it is more convenient if // we swap the leading parts of columns k and r-1 (zero-based) to form // the matrix X_k = [X(0:k, k), X(0:k, r:n)] m_cpqr.m_qr.col(k).head(k + 1).swap( m_cpqr.m_qr.col(rank - 1).head(k + 1)); } // Construct Householder reflector Z(k) to zero out the last row of X_k, // i.e. choose Z(k) such that // [X(k, k), X(k, r:n)] * Z(k) = [beta, 0, .., 0]. RealScalar beta; m_cpqr.m_qr.row(k) .tail(cols - rank + 1) .makeHouseholderInPlace(m_zCoeffs(k), beta); m_cpqr.m_qr(k, rank - 1) = beta; if (k > 0) { // Apply Z(k) to the first k rows of X_k m_cpqr.m_qr.topRightCorner(k, cols - rank + 1) .applyHouseholderOnTheRight( m_cpqr.m_qr.row(k).tail(cols - rank).transpose(), m_zCoeffs(k), &m_temp(0)); } if (k != rank - 1) { // Swap X(0:k,k) back to its proper location. m_cpqr.m_qr.col(k).head(k + 1).swap( m_cpqr.m_qr.col(rank - 1).head(k + 1)); } } } return *this; }
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::dimensionOfKernel | ( | ) | const [inline] |
Definition at line 221 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.dimensionOfKernel(); }
const HCoeffsType& Eigen::CompleteOrthogonalDecomposition< _MatrixType >::hCoeffs | ( | ) | const [inline] |
Q
.For advanced uses only.
Definition at line 268 of file CompleteOrthogonalDecomposition.h.
CompleteOrthogonalDecomposition< MatrixType >::HouseholderSequenceType Eigen::CompleteOrthogonalDecomposition< MatrixType >::householderQ | ( | void | ) | const |
Definition at line 524 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.householderQ(); }
ComputationInfo Eigen::CompleteOrthogonalDecomposition< _MatrixType >::info | ( | ) | const [inline] |
Reports whether the complete orthogonal decomposition was succesful.
Success
. It is provided for compatibility with other factorization routines. Success
Definition at line 342 of file CompleteOrthogonalDecomposition.h.
{ eigen_assert(m_cpqr.m_isInitialized && "Decomposition is not initialized."); return Success; }
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::isInjective | ( | ) | const [inline] |
Definition at line 230 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.isInjective(); }
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::isInvertible | ( | ) | const [inline] |
Definition at line 248 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.isInvertible(); }
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::isSurjective | ( | ) | const [inline] |
Definition at line 239 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.isSurjective(); }
MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< MatrixType >::logAbsDeterminant | ( | ) | const |
Definition at line 375 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.logAbsDeterminant(); }
HouseholderSequenceType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixQ | ( | void | ) | const [inline] |
Definition at line 139 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.householderQ(); }
const MatrixType& Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixQTZ | ( | ) | const [inline] |
Definition at line 152 of file CompleteOrthogonalDecomposition.h.
const MatrixType& Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixT | ( | ) | const [inline] |
matrixT().template triangularView<Upper>()
Definition at line 165 of file CompleteOrthogonalDecomposition.h.
MatrixType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixZ | ( | ) | const [inline] |
Definition at line 143 of file CompleteOrthogonalDecomposition.h.
{ MatrixType Z = MatrixType::Identity(m_cpqr.cols(), m_cpqr.cols()); applyZAdjointOnTheLeftInPlace(Z); return Z.adjoint(); }
RealScalar Eigen::CompleteOrthogonalDecomposition< _MatrixType >::maxPivot | ( | ) | const [inline] |
Definition at line 332 of file CompleteOrthogonalDecomposition.h.
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::nonzeroPivots | ( | ) | const [inline] |
Definition at line 327 of file CompleteOrthogonalDecomposition.h.
{ return m_cpqr.nonzeroPivots(); }
const Inverse<CompleteOrthogonalDecomposition> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::pseudoInverse | ( | ) | const [inline] |
this->pseudoInverse()*rhs
to solve a linear systems. It is more efficient and numerically stable to call this->solve(rhs)
. Definition at line 255 of file CompleteOrthogonalDecomposition.h.
{
return Inverse<CompleteOrthogonalDecomposition>(*this);
}
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::rank | ( | ) | const [inline] |
Definition at line 212 of file CompleteOrthogonalDecomposition.h.
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::rows | ( | void | ) | const [inline] |
Definition at line 260 of file CompleteOrthogonalDecomposition.h.
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.
threshold | The new value to use as the threshold. |
A pivot will be considered nonzero if its absolute value is strictly greater than where maxpivot is the biggest pivot.
If you want to come back to the default behavior, call setThreshold(Default_t)
Definition at line 296 of file CompleteOrthogonalDecomposition.h.
{ m_cpqr.setThreshold(threshold); return *this; }
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 309 of file CompleteOrthogonalDecomposition.h.
{ m_cpqr.setThreshold(Default); return *this; }
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
, where A is the matrix of which *this
is the complete orthogonal decomposition.
B | the right-hand sides of the problem to solve. |
Definition at line 131 of file CompleteOrthogonalDecomposition.h.
{ eigen_assert(m_cpqr.m_isInitialized && "CompleteOrthogonalDecomposition is not initialized."); return Solve<CompleteOrthogonalDecomposition, Rhs>(*this, b.derived()); }
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 318 of file CompleteOrthogonalDecomposition.h.
const HCoeffsType& Eigen::CompleteOrthogonalDecomposition< _MatrixType >::zCoeffs | ( | ) | const [inline] |
Z
.For advanced uses only.
Definition at line 275 of file CompleteOrthogonalDecomposition.h.
{ return m_zCoeffs; }
ColPivHouseholderQR<MatrixType> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::m_cpqr [protected] |
Definition at line 362 of file CompleteOrthogonalDecomposition.h.
RowVectorType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::m_temp [protected] |
Definition at line 364 of file CompleteOrthogonalDecomposition.h.
HCoeffsType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::m_zCoeffs [protected] |
Definition at line 363 of file CompleteOrthogonalDecomposition.h.