Two-sided Jacobi SVD decomposition of a rectangular matrix. More...
#include <JacobiSVD.h>
Two-sided Jacobi SVD decomposition of a rectangular 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:
Definition at line 478 of file JacobiSVD.h.
typedef internal::plain_col_type<MatrixType>::type Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::ColType |
Definition at line 504 of file JacobiSVD.h.
typedef MatrixType::Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::Index |
Definition at line 485 of file JacobiSVD.h.
typedef _MatrixType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::MatrixType |
Definition at line 482 of file JacobiSVD.h.
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::MatrixUType |
Definition at line 498 of file JacobiSVD.h.
typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::MatrixVType |
Definition at line 501 of file JacobiSVD.h.
typedef NumTraits<typename MatrixType::Scalar>::Real Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::RealScalar |
Definition at line 484 of file JacobiSVD.h.
typedef internal::plain_row_type<MatrixType>::type Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::RowType |
Definition at line 503 of file JacobiSVD.h.
typedef MatrixType::Scalar Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::Scalar |
Definition at line 483 of file JacobiSVD.h.
typedef internal::plain_diag_type<MatrixType, RealScalar>::type Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::SingularValuesType |
Definition at line 502 of file JacobiSVD.h.
typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::WorkMatrixType |
Definition at line 507 of file JacobiSVD.h.
anonymous enum |
RowsAtCompileTime | |
ColsAtCompileTime | |
DiagSizeAtCompileTime | |
MaxRowsAtCompileTime | |
MaxColsAtCompileTime | |
MaxDiagSizeAtCompileTime | |
MatrixOptions |
Definition at line 486 of file JacobiSVD.h.
Eigen::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&).
Definition at line 514 of file JacobiSVD.h.
Eigen::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.
Definition at line 528 of file JacobiSVD.h.
Eigen::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.
Definition at line 547 of file JacobiSVD.h.
void Eigen::JacobiSVD< MatrixType, QRPreconditioner >::allocate | ( | Index | rows, |
Index | cols, | ||
unsigned int | computationOptions | ||
) | [private] |
Definition at line 679 of file JacobiSVD.h.
Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::cols | ( | void | ) | const [inline] |
Definition at line 653 of file JacobiSVD.h.
JacobiSVD< MatrixType, QRPreconditioner > & Eigen::JacobiSVD< MatrixType, QRPreconditioner >::compute | ( | const MatrixType & | matrix, |
unsigned int | computationOptions | ||
) |
Method performing the decomposition of given matrix using custom options.
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.
Definition at line 726 of file JacobiSVD.h.
JacobiSVD& Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::compute | ( | const MatrixType & | matrix | ) | [inline] |
Method performing the decomposition of given matrix using current options.
matrix | the matrix to decompose |
This method uses the current computationOptions, as already passed to the constructor or to compute(const MatrixType&, unsigned int).
Definition at line 574 of file JacobiSVD.h.
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::computeU | ( | ) | const [inline] |
Definition at line 623 of file JacobiSVD.h.
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::computeV | ( | ) | const [inline] |
Definition at line 625 of file JacobiSVD.h.
const MatrixUType& Eigen::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.
Definition at line 588 of file JacobiSVD.h.
const MatrixVType& Eigen::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.
Definition at line 604 of file JacobiSVD.h.
Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::nonzeroSingularValues | ( | ) | const [inline] |
Definition at line 646 of file JacobiSVD.h.
Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::rows | ( | void | ) | const [inline] |
Definition at line 652 of file JacobiSVD.h.
const SingularValuesType& Eigen::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.
Definition at line 616 of file JacobiSVD.h.
const internal::solve_retval<JacobiSVD, Rhs> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::solve | ( | const MatrixBase< Rhs > & | b | ) | const [inline] |
b | the right-hand-side of the equation to solve. |
Definition at line 638 of file JacobiSVD.h.
friend struct internal::qr_preconditioner_impl [friend] |
Definition at line 672 of file JacobiSVD.h.
friend struct internal::svd_precondition_2x2_block_to_be_real [friend] |
Definition at line 670 of file JacobiSVD.h.
Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_cols [protected] |
Definition at line 667 of file JacobiSVD.h.
unsigned int Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_computationOptions [protected] |
Definition at line 666 of file JacobiSVD.h.
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_computeFullU [protected] |
Definition at line 664 of file JacobiSVD.h.
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_computeFullV [protected] |
Definition at line 665 of file JacobiSVD.h.
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_computeThinU [protected] |
Definition at line 664 of file JacobiSVD.h.
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_computeThinV [protected] |
Definition at line 665 of file JacobiSVD.h.
Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_diagSize [protected] |
Definition at line 667 of file JacobiSVD.h.
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_isAllocated [protected] |
Definition at line 663 of file JacobiSVD.h.
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_isInitialized [protected] |
Definition at line 663 of file JacobiSVD.h.
MatrixUType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_matrixU [protected] |
Definition at line 659 of file JacobiSVD.h.
MatrixVType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_matrixV [protected] |
Definition at line 660 of file JacobiSVD.h.
Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_nonzeroSingularValues [protected] |
Definition at line 667 of file JacobiSVD.h.
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_qr_precond_morecols [protected] |
Definition at line 674 of file JacobiSVD.h.
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_qr_precond_morerows [protected] |
Definition at line 675 of file JacobiSVD.h.
Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_rows [protected] |
Definition at line 667 of file JacobiSVD.h.
SingularValuesType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_singularValues [protected] |
Definition at line 661 of file JacobiSVD.h.
WorkMatrixType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_workMatrix [protected] |
Definition at line 662 of file JacobiSVD.h.