$darkmode
#include <Eigen/src/SVD/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 |
| Options | this optional parameter allows one to specify the type of QR decomposition that will be used internally for the R-SVD step for non-square matrices. Additionally, it allows one to specify whether to compute thin or full unitaries U and V. See discussion of possible values below. |
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.
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:
Here is the matrix m: -0.824 -0.199 0.924 -0.237 -0.0532 0.146 Its singular values are: 1.24 0.338 Its left singular vectors are the columns of the thin U matrix: -0.657 -0.695 0.753 -0.579 -0.048 0.426 Its right singular vectors are the columns of the thin V matrix: 0.999 0.044 -0.044 0.999 Now consider this rhs vector: 1 0 0 A least-squares solution of m*x = rhs is: -0.62 -2.03
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 \( O(n^2p) \) 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 QR preconditioners that can be set with Options template parameter are:
One may also use the Options template parameter to specify how the unitaries should be computed. The options are ComputeThinU, ComputeThinV, ComputeFullU, ComputeFullV. It is not possible to request both the thin and full versions of a unitary. By default, unitaries will not be computed.
You can set the QRPreconditioner and unitary options together: JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner | ComputeThinU | ComputeFullV>
Inheritance diagram for Eigen::JacobiSVD< MatrixType_, Options_ >:
Public Member Functions | |
| template<typename Derived > | |
| JacobiSVD & | compute (const MatrixBase< Derived > &matrix) |
| Method performing the decomposition of given matrix. Computes Thin/Full unitaries U/V if specified using the Options template parameter or the class constructor. More... | |
| template<typename Derived > | |
| JacobiSVD & | compute (const MatrixBase< Derived > &matrix, unsigned int computationOptions) |
Method performing the decomposition of given matrix, as specified by the computationOptions parameter. More... | |
| bool | computeU () const |
| bool | computeV () const |
| JacobiSVD () | |
| Default Constructor. More... | |
| JacobiSVD (Index rows, Index cols) | |
| Default Constructor with memory preallocation. More... | |
| JacobiSVD (Index rows, Index cols, unsigned int computationOptions) | |
| Default Constructor with memory preallocation. More... | |
| template<typename Derived > | |
| JacobiSVD (const MatrixBase< Derived > &matrix) | |
| Constructor performing the decomposition of given matrix, using the custom options specified with the Options template parameter. More... | |
| template<typename Derived > | |
| JacobiSVD (const MatrixBase< Derived > &matrix, unsigned int computationOptions) | |
| Constructor performing the decomposition of given matrix using specified options for computing unitaries. More... | |
| Index | rank () const |
Public Member Functions inherited from Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > > | |
| bool | computeU () const |
| bool | computeV () const |
| ComputationInfo | info () const |
| Reports whether previous computation was successful. More... | |
| const MatrixUType & | matrixU () const |
| const MatrixVType & | matrixV () const |
| Index | nonzeroSingularValues () const |
| Index | rank () const |
| JacobiSVD< MatrixType_, Options_ > & | setThreshold (const RealScalar &threshold) |
| JacobiSVD< MatrixType_, Options_ > & | setThreshold (Default_t) |
| const SingularValuesType & | singularValues () const |
| const Solve< JacobiSVD< MatrixType_, Options_ >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| RealScalar | threshold () const |
Public Member Functions inherited from Eigen::SolverBase< SVDBase< JacobiSVD< MatrixType_, Options_ > > > | |
| const AdjointReturnType | adjoint () const |
| constexpr const SVDBase< JacobiSVD< MatrixType_, Options_ > > & | derived () const |
| constexpr SVDBase< JacobiSVD< MatrixType_, Options_ > > & | derived () |
| const Solve< SVDBase< JacobiSVD< MatrixType_, Options_ > >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| SolverBase () | |
| const ConstTransposeReturnType | transpose () const |
Public Member Functions inherited from Eigen::EigenBase< Derived > | |
| constexpr Index | cols () const noexcept |
| constexpr Derived & | derived () |
| constexpr const Derived & | derived () const |
| constexpr Index | rows () const noexcept |
| constexpr Index | size () const noexcept |
Additional Inherited Members | |
Public Types inherited from Eigen::EigenBase< Derived > | |
| typedef Eigen::Index | Index |
| The interface type of indices. More... | |
Protected Member Functions inherited from Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > > | |
| SVDBase () | |
| Default Constructor. More... | |
|
inline |
Default Constructor.
The default constructor is useful in cases in which the user intends to perform decompositions via JacobiSVD::compute(const MatrixType&).
|
inline |
Default Constructor with memory preallocation.
Like the default constructor but with preallocation of the internal data according to the specified problem size and Options template parameter.
|
inline |
Default Constructor with memory preallocation.
Like the default constructor but with preallocation of the internal data according to the specified problem size.
One cannot request unitaries using both the Options template parameter and the constructor. If possible, prefer using the Options template parameter.
| rows | number of rows for the input matrix |
| cols | number of columns for the input matrix |
| computationOptions | specify whether to compute Thin/Full unitaries U/V |
|
inlineexplicit |
Constructor performing the decomposition of given matrix, using the custom options specified with the Options template parameter.
| matrix | the matrix to decompose |
|
inline |
Constructor performing the decomposition of given matrix using specified options for computing unitaries.
One cannot request unitiaries using both the Options template parameter and the constructor. If possible, prefer using the Options template parameter.
| matrix | the matrix to decompose |
| computationOptions | specify whether to compute Thin/Full unitaries U/V |
|
inline |
Method performing the decomposition of given matrix. Computes Thin/Full unitaries U/V if specified using the Options template parameter or the class constructor.
| matrix | the matrix to decompose |
|
inline |
Method performing the decomposition of given matrix, as specified by the computationOptions parameter.
| matrix | the matrix to decompose |
| computationOptions | specify whether to compute Thin/Full unitaries U/V |
|
inline |
|
inline |
|
inline |
*this is the SVD.