Eigen  3.3.3
Eigen::SVDBase< Derived > Class Template Reference

Detailed Description

template<typename Derived>
class Eigen::SVDBase< Derived >

Base class of SVD algorithms.

Template Parameters:
Derivedthe type of the actual SVD decomposition

SVD decomposition consists in decomposing any n-by-p matrix A as a product

\[ A = U S V^* \]

where U is a n-by-n unitary, V is a p-by-p unitary, and S is a n-by-p real positive matrix which is zero outside of its main diagonal; the diagonal entries of S are known as the singular values of A and the columns of U and V are known as the left and right singular vectors of A respectively.

Singular values are always sorted in decreasing order.

You can ask for only thin U or V to be computed, meaning the following. In case of a rectangular n-by-p matrix, letting m be the smaller value among n and p, there are only m singular vectors; the remaining columns of U and V do not correspond to actual singular vectors. Asking for thin U or V means asking for only their m first columns to be formed. So U is then a n-by-m matrix, and V is then a p-by-m matrix. Notice that thin U and V are all you need for (least squares) solving.

If the input matrix has inf or nan coefficients, the result of the computation is undefined, but the computation is guaranteed to terminate in finite (and reasonable) time.

See also:
class BDCSVD, class JacobiSVD

List of all members.

Public Types

typedef Eigen::Index Index

Public Member Functions

bool computeU () const
bool computeV () const
const MatrixUTypematrixU () const
const MatrixVTypematrixV () const
Index nonzeroSingularValues () const
Index rank () const
Derived & setThreshold (const RealScalar &threshold)
Derived & setThreshold (Default_t)
const SingularValuesType & singularValues () const
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
RealScalar threshold () const

Protected Member Functions

 SVDBase ()
 Default Constructor.

Member Typedef Documentation

template<typename Derived>
typedef Eigen::Index Eigen::SVDBase< Derived >::Index
Deprecated:
since Eigen 3.3

Constructor & Destructor Documentation

template<typename Derived>
Eigen::SVDBase< Derived >::SVDBase ( ) [inline, protected]

Default Constructor.

Default constructor of SVDBase


Member Function Documentation

template<typename Derived>
bool Eigen::SVDBase< Derived >::computeU ( ) const [inline]
Returns:
true if U (full or thin) is asked for in this SVD decomposition
template<typename Derived>
bool Eigen::SVDBase< Derived >::computeV ( ) const [inline]
Returns:
true if V (full or thin) is asked for in this SVD decomposition
template<typename Derived>
const MatrixUType& Eigen::SVDBase< Derived >::matrixU ( ) const [inline]
Returns:
the U matrix.

For the SVD decomposition of a n-by-p matrix, letting m be the minimum of n and p, the U matrix is n-by-n if you asked for ComputeFullU , and is n-by-m if you asked for ComputeThinU .

The m first columns of U are the left singular vectors of the matrix being decomposed.

This method asserts that you asked for U to be computed.

template<typename Derived>
const MatrixVType& Eigen::SVDBase< Derived >::matrixV ( ) const [inline]
Returns:
the V matrix.

For the SVD decomposition of a n-by-p matrix, letting m be the minimum of n and p, the V matrix is p-by-p if you asked for ComputeFullV , and is p-by-m if you asked for ComputeThinV .

The m first columns of V are the right singular vectors of the matrix being decomposed.

This method asserts that you asked for V to be computed.

template<typename Derived>
Index Eigen::SVDBase< Derived >::nonzeroSingularValues ( ) const [inline]
Returns:
the number of singular values that are not exactly 0
template<typename Derived>
Index Eigen::SVDBase< Derived >::rank ( ) const [inline]
Returns:
the rank of the matrix of which *this is the SVD.
Note:
This method has to determine which singular values should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).
template<typename Derived>
Derived& Eigen::SVDBase< Derived >::setThreshold ( const RealScalar &  threshold) [inline]

Allows to prescribe a threshold to be used by certain methods, such as rank() and solve(), which need to determine when singular values are to be considered nonzero. This is not used for the SVD decomposition itself.

When it needs to get the threshold value, Eigen calls threshold(). The default is NumTraits<Scalar>::epsilon()

Parameters:
thresholdThe new value to use as the threshold.

A singular value will be considered nonzero if its value is strictly greater than $ \vert singular value \vert \leqslant threshold \times \vert max singular value \vert $.

If you want to come back to the default behavior, call setThreshold(Default_t)

template<typename Derived>
Derived& Eigen::SVDBase< Derived >::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.

 svd.setThreshold(Eigen::Default); 

See the documentation of setThreshold(const RealScalar&).

template<typename Derived>
const SingularValuesType& Eigen::SVDBase< Derived >::singularValues ( ) const [inline]
Returns:
the vector of singular values.

For the SVD decomposition of a n-by-p matrix, letting m be the minimum of n and p, the returned vector has size m. Singular values are always sorted in decreasing order.

template<typename Derived>
template<typename Rhs >
const Solve<Derived, Rhs> Eigen::SVDBase< Derived >::solve ( const MatrixBase< Rhs > &  b) const [inline]
Returns:
a (least squares) solution of $ A x = b $ using the current SVD decomposition of A.
Parameters:
bthe right-hand-side of the equation to solve.
Note:
Solving requires both U and V to be computed. Thin U and V are enough, there is no need for full U or V.
SVD solving is implicitly least-squares. Thus, this method serves both purposes of exact solving and least-squares solving. In other words, the returned solution is guaranteed to minimize the Euclidean norm $ \Vert A x - b \Vert $.
template<typename Derived>
RealScalar Eigen::SVDBase< Derived >::threshold ( ) const [inline]

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

See the documentation of setThreshold(const RealScalar&).


The documentation for this class was generated from the following file:
 All Classes Functions Variables Typedefs Enumerations Enumerator Friends