Base class of SVD algorithms. More...
#include <SVDBase.h>
Public Types | |
enum | { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime), MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime, MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime), MatrixOptions = MatrixType::Options } |
typedef Eigen::Index | Index |
typedef internal::traits< Derived >::MatrixType | MatrixType |
typedef Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime > | MatrixUType |
typedef Matrix< Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime > | MatrixVType |
typedef NumTraits< typename MatrixType::Scalar >::Real | RealScalar |
typedef MatrixType::Scalar | Scalar |
typedef internal::plain_diag_type< MatrixType, RealScalar >::type | SingularValuesType |
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 |
Index | cols () const |
bool | computeU () const |
bool | computeV () const |
Derived & | derived () |
const Derived & | derived () const |
const MatrixUType & | matrixU () const |
const MatrixVType & | matrixV () const |
Index | nonzeroSingularValues () const |
Index | rank () const |
Index | rows () 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 | |
bool | allocate (Index rows, Index cols, unsigned int computationOptions) |
SVDBase () | |
Default Constructor. More... | |
Static Protected Member Functions | |
static void | check_template_parameters () |
Protected Attributes | |
Index | m_cols |
unsigned int | m_computationOptions |
bool | m_computeFullU |
bool | m_computeFullV |
bool | m_computeThinU |
bool | m_computeThinV |
Index | m_diagSize |
bool | m_isAllocated |
bool | m_isInitialized |
MatrixUType | m_matrixU |
MatrixVType | m_matrixV |
Index | m_nonzeroSingularValues |
RealScalar | m_prescribedThreshold |
Index | m_rows |
SingularValuesType | m_singularValues |
bool | m_usePrescribedThreshold |
Base class of SVD algorithms.
Derived | the type of the actual SVD decomposition |
SVD decomposition consists in decomposing any n-by-p matrix A as a product
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.
typedef Eigen::Index Eigen::SVDBase< Derived >::Index |
typedef internal::traits<Derived>::MatrixType Eigen::SVDBase< Derived >::MatrixType |
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime> Eigen::SVDBase< Derived >::MatrixUType |
typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime> Eigen::SVDBase< Derived >::MatrixVType |
typedef NumTraits<typename MatrixType::Scalar>::Real Eigen::SVDBase< Derived >::RealScalar |
typedef MatrixType::Scalar Eigen::SVDBase< Derived >::Scalar |
typedef internal::plain_diag_type<MatrixType, RealScalar>::type Eigen::SVDBase< Derived >::SingularValuesType |
typedef MatrixType::StorageIndex Eigen::SVDBase< Derived >::StorageIndex |
anonymous enum |
|
inlineprotected |
EIGEN_DEVICE_FUNC void Eigen::SVDBase< Derived >::_solve_impl | ( | const RhsType & | rhs, |
DstType & | dst | ||
) | const |
void Eigen::SVDBase< Derived >::_solve_impl | ( | const RhsType & | rhs, |
DstType & | dst | ||
) | const |
|
protected |
|
inlinestaticprotected |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
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.
|
inline |
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.
|
inline |
|
inline |
*this
is the SVD.
|
inline |
|
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()
threshold | The new value to use as the threshold. |
A singular value will be considered nonzero if its value is strictly greater than .
If you want to come back to the default behavior, call 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.
See the documentation of setThreshold(const RealScalar&).
|
inline |
|
inline |
b | the right-hand-side of the equation to solve. |
|
inline |
Returns the threshold that will be used by certain methods such as rank().
See the documentation of setThreshold(const RealScalar&).
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |