Two-sided Jacobi SVD decomposition of a square matrix. More...
#include <JacobiSVD.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 _MatrixType | MatrixType |
typedef MatrixType::Scalar | Scalar |
typedef NumTraits< typename MatrixType::Scalar >::Real | RealScalar |
typedef MatrixType::Index | Index |
typedef Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime > | MatrixUType |
typedef Matrix< Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime > | MatrixVType |
typedef internal::plain_diag_type < MatrixType, RealScalar > ::type | SingularValuesType |
typedef internal::plain_row_type < MatrixType >::type | RowType |
typedef internal::plain_col_type < MatrixType >::type | ColType |
typedef Matrix< Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime > | WorkMatrixType |
Public Member Functions | |
JacobiSVD () | |
Default Constructor. | |
JacobiSVD (Index rows, Index cols, unsigned int computationOptions=0) | |
Default Constructor with memory preallocation. | |
JacobiSVD (const MatrixType &matrix, unsigned int computationOptions=0) | |
Constructor performing the decomposition of given matrix. | |
JacobiSVD & | compute (const MatrixType &matrix, unsigned int computationOptions=0) |
Method performing the decomposition of given matrix. | |
const MatrixUType & | matrixU () const |
const MatrixVType & | matrixV () const |
const SingularValuesType & | singularValues () const |
bool | computeU () const |
bool | computeV () const |
template<typename Rhs > | |
const internal::solve_retval < JacobiSVD, Rhs > | solve (const MatrixBase< Rhs > &b) const |
Index | nonzeroSingularValues () const |
Index | rows () const |
Index | cols () const |
Protected Attributes | |
MatrixUType | m_matrixU |
MatrixVType | m_matrixV |
SingularValuesType | m_singularValues |
WorkMatrixType | m_workMatrix |
bool | m_isInitialized |
bool | m_computeFullU |
bool | m_computeThinU |
bool | m_computeFullV |
bool | m_computeThinV |
Index | m_nonzeroSingularValues |
Index | m_rows |
Index | m_cols |
Index | m_diagSize |
Friends | |
struct | internal::svd_precondition_2x2_block_to_be_real |
struct | internal::qr_preconditioner_impl |
Two-sided Jacobi SVD decomposition of a square matrix.
MatrixType | the type of the matrix of which we are computing the SVD decomposition | |
QRPreconditioner | this optional parameter allows to specify the type of QR decomposition that will be used internally for the R-SVD step for non-square matrices. See discussion of possible values below. |
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.
This JacobiSVD decomposition computes only the singular values by default. If you want U or V, you need to ask for them explicitly.
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.
Here's an example demonstrating basic usage:
Output:
This JacobiSVD class is a two-sided Jacobi R-SVD decomposition, ensuring optimal reliability and accuracy. The downside is that it's slower than bidiagonalizing SVD algorithms for large square matrices; however its complexity is still where n is the smaller dimension and p is the greater dimension, meaning that it is still of the same order of complexity as the faster bidiagonalizing R-SVD algorithms. In particular, like any R-SVD, it takes advantage of non-squareness in that its complexity is only linear in the greater dimension.
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.
The possible values for QRPreconditioner are:
typedef internal::plain_col_type<MatrixType>::type JacobiSVD< _MatrixType, QRPreconditioner >::ColType |
typedef MatrixType::Index JacobiSVD< _MatrixType, QRPreconditioner >::Index |
typedef _MatrixType JacobiSVD< _MatrixType, QRPreconditioner >::MatrixType |
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime> JacobiSVD< _MatrixType, QRPreconditioner >::MatrixUType |
typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime> JacobiSVD< _MatrixType, QRPreconditioner >::MatrixVType |
typedef NumTraits<typename MatrixType::Scalar>::Real JacobiSVD< _MatrixType, QRPreconditioner >::RealScalar |
typedef internal::plain_row_type<MatrixType>::type JacobiSVD< _MatrixType, QRPreconditioner >::RowType |
typedef MatrixType::Scalar JacobiSVD< _MatrixType, QRPreconditioner >::Scalar |
typedef internal::plain_diag_type<MatrixType, RealScalar>::type JacobiSVD< _MatrixType, QRPreconditioner >::SingularValuesType |
typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime> JacobiSVD< _MatrixType, QRPreconditioner >::WorkMatrixType |
anonymous enum |
JacobiSVD< _MatrixType, QRPreconditioner >::JacobiSVD | ( | ) | [inline] |
Default Constructor.
The default constructor is useful in cases in which the user intends to perform decompositions via JacobiSVD::compute(const MatrixType&).
JacobiSVD< _MatrixType, QRPreconditioner >::JacobiSVD | ( | Index | rows, | |
Index | cols, | |||
unsigned int | computationOptions = 0 | |||
) | [inline] |
Default Constructor with memory preallocation.
Like the default constructor but with preallocation of the internal data according to the specified problem size.
JacobiSVD< _MatrixType, QRPreconditioner >::JacobiSVD | ( | const MatrixType & | matrix, | |
unsigned int | computationOptions = 0 | |||
) | [inline] |
Constructor performing the decomposition of given matrix.
matrix | the matrix to decompose | |
computationOptions | optional parameter allowing to specify if you want full or thin U or V unitaries to be computed. By default, none is computed. This is a bit-field, the possible bits are ComputeFullU, ComputeThinU, ComputeFullV, ComputeThinV. |
Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not available with the (non-default) FullPivHouseholderQR preconditioner.
Index JacobiSVD< _MatrixType, QRPreconditioner >::cols | ( | void | ) | const [inline] |
JacobiSVD< MatrixType, QRPreconditioner > & JacobiSVD< MatrixType, QRPreconditioner >::compute | ( | const MatrixType & | matrix, | |
unsigned int | computationOptions = 0 | |||
) |
Method performing the decomposition of given matrix.
matrix | the matrix to decompose | |
computationOptions | optional parameter allowing to specify if you want full or thin U or V unitaries to be computed. By default, none is computed. This is a bit-field, the possible bits are ComputeFullU, ComputeThinU, ComputeFullV, ComputeThinV. |
Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not available with the (non-default) FullPivHouseholderQR preconditioner.
bool JacobiSVD< _MatrixType, QRPreconditioner >::computeU | ( | ) | const [inline] |
bool JacobiSVD< _MatrixType, QRPreconditioner >::computeV | ( | ) | const [inline] |
const MatrixUType& JacobiSVD< _MatrixType, QRPreconditioner >::matrixU | ( | ) | const [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.
const MatrixVType& JacobiSVD< _MatrixType, QRPreconditioner >::matrixV | ( | ) | const [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.
Index JacobiSVD< _MatrixType, QRPreconditioner >::nonzeroSingularValues | ( | ) | const [inline] |
Index JacobiSVD< _MatrixType, QRPreconditioner >::rows | ( | void | ) | const [inline] |
const SingularValuesType& JacobiSVD< _MatrixType, QRPreconditioner >::singularValues | ( | ) | const [inline] |
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.
const internal::solve_retval<JacobiSVD, Rhs> JacobiSVD< _MatrixType, QRPreconditioner >::solve | ( | const MatrixBase< Rhs > & | b | ) | const [inline] |
b | the right-hand-side of the equation to solve. |
friend struct internal::qr_preconditioner_impl [friend] |
friend struct internal::svd_precondition_2x2_block_to_be_real [friend] |
Index JacobiSVD< _MatrixType, QRPreconditioner >::m_cols [protected] |
bool JacobiSVD< _MatrixType, QRPreconditioner >::m_computeFullU [protected] |
bool JacobiSVD< _MatrixType, QRPreconditioner >::m_computeFullV [protected] |
bool JacobiSVD< _MatrixType, QRPreconditioner >::m_computeThinU [protected] |
bool JacobiSVD< _MatrixType, QRPreconditioner >::m_computeThinV [protected] |
Index JacobiSVD< _MatrixType, QRPreconditioner >::m_diagSize [protected] |
bool JacobiSVD< _MatrixType, QRPreconditioner >::m_isInitialized [protected] |
MatrixUType JacobiSVD< _MatrixType, QRPreconditioner >::m_matrixU [protected] |
MatrixVType JacobiSVD< _MatrixType, QRPreconditioner >::m_matrixV [protected] |
Index JacobiSVD< _MatrixType, QRPreconditioner >::m_nonzeroSingularValues [protected] |
Index JacobiSVD< _MatrixType, QRPreconditioner >::m_rows [protected] |
SingularValuesType JacobiSVD< _MatrixType, QRPreconditioner >::m_singularValues [protected] |
WorkMatrixType JacobiSVD< _MatrixType, QRPreconditioner >::m_workMatrix [protected] |