Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
Eigen::JacobiSVD< MatrixType_, Options_ > Class Template Reference

Detailed Description

template<typename MatrixType_, int Options_>
class Eigen::JacobiSVD< MatrixType_, Options_ >

Two-sided Jacobi SVD decomposition of a rectangular matrix.

Template Parameters
MatrixType_the type of the matrix of which we are computing the SVD decomposition
Optionsthis 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:

cout << "Here is the matrix m:" << endl << m << endl;
JacobiSVD<MatrixXf, ComputeThinU | ComputeThinV> svd(m);
cout << "Its singular values are:" << endl << svd.singularValues() << endl;
cout << "Its left singular vectors are the columns of the thin U matrix:" << endl << svd.matrixU() << endl;
cout << "Its right singular vectors are the columns of the thin V matrix:" << endl << svd.matrixV() << endl;
Vector3f rhs(1, 0, 0);
cout << "Now consider this rhs vector:" << endl << rhs << endl;
cout << "A least-squares solution of m*x = rhs is:" << endl << svd.solve(rhs) << endl;
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< float, 3, 1 > Vector3f
3×1 vector of type float.
Definition: Matrix.h:500
Matrix< float, Dynamic, Dynamic > MatrixXf
Dynamic×Dynamic matrix of type float.
Definition: Matrix.h:500

Output:

Here is the matrix m:
  0.68  0.597
-0.211  0.823
 0.566 -0.605
Its singular values are:
 1.19
0.899
Its left singular vectors are the columns of the thin U matrix:
  0.388   0.866
  0.712 -0.0634
 -0.586   0.496
Its right singular vectors are the columns of the thin V matrix:
-0.183  0.983
 0.983  0.183
Now consider this rhs vector:
1
0
0
A least-squares solution of m*x = rhs is:
0.888
0.496

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:

  • ColPivHouseholderQRPreconditioner is the default. In practice it's very safe. It uses column-pivoting QR.
  • FullPivHouseholderQRPreconditioner, is the safest and slowest. It uses full-pivoting QR. Contrary to other QRs, it doesn't allow computing thin unitaries.
  • HouseholderQRPreconditioner is the fastest, and less safe and accurate than the pivoting variants. It uses non-pivoting QR. This is very similar in safety and accuracy to the bidiagonalization process used by bidiagonalizing SVD algorithms (since bidiagonalization is inherently non-pivoting). However the resulting SVD is still more reliable than bidiagonalizing SVDs because the Jacobi-based iterarive process is more reliable than the optimized bidiagonal SVD iterations.
  • NoQRPreconditioner allows not to use a QR preconditioner at all. This is useful if you know that you will only be computing JacobiSVD decompositions of square matrices. Non-square matrices require a QR preconditioner. Using this option will result in faster compilation and smaller executable code. It won't significantly speed up computation, since JacobiSVD is always checking if QR preconditioning is needed before applying it anyway.

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>

See also
MatrixBase::jacobiSvd()
+ Inheritance diagram for Eigen::JacobiSVD< MatrixType_, Options_ >:

Public Member Functions

EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
JacobiSVDcompute (const MatrixType &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...
 
EIGEN_DEPRECATED JacobiSVDcompute (const MatrixType &matrix, unsigned int computationOptions)
 Method performing the decomposition of given matrix, as specified by the computationOptions parameter. More...
 
 JacobiSVD ()
 Default Constructor. More...
 
 JacobiSVD (const MatrixType &matrix)
 Constructor performing the decomposition of given matrix, using the custom options specified with the Options template paramter. More...
 
 JacobiSVD (const MatrixType &matrix, unsigned int computationOptions)
 Constructor performing the decomposition of given matrix using specified options for computing unitaries. More...
 
 JacobiSVD (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
EIGEN_DEPRECATED JacobiSVD (Index rows, Index cols, unsigned int computationOptions)
 Default Constructor with memory preallocation. More...
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
- 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 MatrixUTypematrixU () const
 
const MatrixVTypematrixV () 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< Derived >
const AdjointReturnType adjoint () const
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & derived ()
 
const Derived & derived () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 

Additional Inherited Members

- Public Types inherited from Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >
typedef Eigen::Index Index
 
- 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...
 

Constructor & Destructor Documentation

◆ JacobiSVD() [1/5]

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::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() [2/5]

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD ( Index  rows,
Index  cols 
)
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.

See also
JacobiSVD()

◆ JacobiSVD() [3/5]

template<typename MatrixType_ , int Options_>
EIGEN_DEPRECATED Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD ( Index  rows,
Index  cols,
unsigned int  computationOptions 
)
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.

Parameters
computationOptionsspecify whether to compute Thin/Full unitaries U/V
See also
JacobiSVD()
Deprecated:
Will be removed in the next major Eigen version. Options should be specified in the Options template parameter.

◆ JacobiSVD() [4/5]

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD ( const MatrixType &  matrix)
inlineexplicit

Constructor performing the decomposition of given matrix, using the custom options specified with the Options template paramter.

Parameters
matrixthe matrix to decompose

◆ JacobiSVD() [5/5]

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD ( const MatrixType &  matrix,
unsigned int  computationOptions 
)
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.

Parameters
matrixthe matrix to decompose
computationOptionsspecify whether to compute Thin/Full unitaries U/V
Deprecated:
Will be removed in the next major Eigen version. Options should be specified in the Options template parameter.

Member Function Documentation

◆ cols()

template<typename MatrixType_ , int Options_>
EIGEN_CONSTEXPR Index Eigen::EigenBase< Derived >::cols ( void  )
inline
Returns
the number of columns.
See also
rows(), ColsAtCompileTime

◆ compute() [1/2]

template<typename MatrixType_ , int Options_>
JacobiSVD& Eigen::JacobiSVD< MatrixType_, Options_ >::compute ( const MatrixType &  matrix)
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.

Parameters
matrixthe matrix to decompose

◆ compute() [2/2]

template<typename MatrixType_ , int Options_>
EIGEN_DEPRECATED JacobiSVD& Eigen::JacobiSVD< MatrixType_, Options_ >::compute ( const MatrixType &  matrix,
unsigned int  computationOptions 
)
inline

Method performing the decomposition of given matrix, as specified by the computationOptions parameter.

Parameters
matrixthe matrix to decompose
computationOptionsspecify whether to compute Thin/Full unitaries U/V
Deprecated:
Will be removed in the next major Eigen version. Options should be specified in the Options template parameter.

◆ rows()

template<typename MatrixType_ , int Options_>
EIGEN_CONSTEXPR Index Eigen::EigenBase< Derived >::rows ( void  )
inline
Returns
the number of rows.
See also
cols(), RowsAtCompileTime

The documentation for this class was generated from the following file: