template<typename _MatrixType, int QRPreconditioner>
Eigen::JacobiSVD class

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
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

\[ 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:

MatrixXf m = MatrixXf::Random(3,2);
cout << "Here is the matrix m:" << endl << m << endl;
JacobiSVD<MatrixXf> svd(m, ComputeThinU | ComputeThinV);
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;

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 values for QRPreconditioner 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.

Base classes

template<typename Derived>
class SVDBase
Base class of SVD algorithms.

Public types

enum (anonymous) { 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 }
using ColType = internal::plain_col_type<MatrixType>::type
using MatrixType = _MatrixType
using MatrixUType = Base::MatrixUType
using MatrixVType = Base::MatrixVType
using RealScalar = NumTraits<typename MatrixType::Scalar>::Real
using RowType = internal::plain_row_type<MatrixType>::type
using Scalar = MatrixType::Scalar
using SingularValuesType = Base::SingularValuesType
using WorkMatrixType = Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>

Constructors, destructors, conversion operators

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) explicit
Constructor performing the decomposition of given matrix.

Public functions

auto cols() const -> Index
auto compute(const MatrixType& matrix, unsigned int computationOptions) -> JacobiSVD&
Method performing the decomposition of given matrix using custom options.
auto compute(const MatrixType& matrix) -> JacobiSVD&
Method performing the decomposition of given matrix using current options.
auto computeU() const -> bool
auto computeV() const -> bool
auto rank() const -> Index
auto rows() const -> Index

Protected variables

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
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows
Index m_rows
MatrixType m_scaledMatrix
SingularValuesType m_singularValues
bool m_usePrescribedThreshold
WorkMatrixType m_workMatrix

Enum documentation

template<typename _MatrixType, int QRPreconditioner>
enum Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::(anonymous)

Enumerators
RowsAtCompileTime
ColsAtCompileTime
DiagSizeAtCompileTime
MaxRowsAtCompileTime
MaxColsAtCompileTime
MaxDiagSizeAtCompileTime
MatrixOptions

Typedef documentation

template<typename _MatrixType, int QRPreconditioner>
typedef internal::plain_col_type<MatrixType>::type Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::ColType

template<typename _MatrixType, int QRPreconditioner>
typedef _MatrixType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::MatrixType

template<typename _MatrixType, int QRPreconditioner>
typedef Base::MatrixUType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::MatrixUType

template<typename _MatrixType, int QRPreconditioner>
typedef Base::MatrixVType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::MatrixVType

template<typename _MatrixType, int QRPreconditioner>
typedef NumTraits<typename MatrixType::Scalar>::Real Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::RealScalar

template<typename _MatrixType, int QRPreconditioner>
typedef internal::plain_row_type<MatrixType>::type Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::RowType

template<typename _MatrixType, int QRPreconditioner>
typedef MatrixType::Scalar Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::Scalar

template<typename _MatrixType, int QRPreconditioner>
typedef Base::SingularValuesType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::SingularValuesType

template<typename _MatrixType, int QRPreconditioner>
typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime> Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::WorkMatrixType

Function documentation

template<typename _MatrixType, int QRPreconditioner>
Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::JacobiSVD()

Default Constructor.

The default constructor is useful in cases in which the user intends to perform decompositions via JacobiSVD::compute(const MatrixType&).

template<typename _MatrixType, int QRPreconditioner>
Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::JacobiSVD(Index rows, Index cols, unsigned int computationOptions = 0)

Default Constructor with memory preallocation.

Like the default constructor but with preallocation of the internal data according to the specified problem size.

template<typename _MatrixType, int QRPreconditioner>
Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::JacobiSVD(const MatrixType& matrix, unsigned int computationOptions = 0) explicit

Constructor performing the decomposition of given matrix.

Parameters
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.

template<typename _MatrixType, int QRPreconditioner>
Index Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::cols() const

template<typename _MatrixType, int QRPreconditioner>
JacobiSVD& Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)

Method performing the decomposition of given matrix using custom options.

Parameters
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.

template<typename _MatrixType, int QRPreconditioner>
JacobiSVD& Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::compute(const MatrixType& matrix)

Method performing the decomposition of given matrix using current options.

Parameters
matrix the matrix to decompose

This method uses the current computationOptions, as already passed to the constructor or to compute(const MatrixType&, unsigned int).

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::computeU() const

Returns true if U (full or thin) is asked for in this SVD decomposition

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::computeV() const

Returns true if V (full or thin) is asked for in this SVD decomposition

template<typename _MatrixType, int QRPreconditioner>
Index Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::rank() const

Returns the rank of the matrix of which *this is the SVD.

template<typename _MatrixType, int QRPreconditioner>
Index Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::rows() const

Variable documentation

template<typename _MatrixType, int QRPreconditioner>
Index Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_cols protected

template<typename _MatrixType, int QRPreconditioner>
unsigned int Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_computationOptions protected

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_computeFullU protected

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_computeFullV protected

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_computeThinU protected

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_computeThinV protected

template<typename _MatrixType, int QRPreconditioner>
Index Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_diagSize protected

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_isAllocated protected

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_isInitialized protected

template<typename _MatrixType, int QRPreconditioner>
MatrixUType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_matrixU protected

template<typename _MatrixType, int QRPreconditioner>
MatrixVType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_matrixV protected

template<typename _MatrixType, int QRPreconditioner>
Index Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_nonzeroSingularValues protected

template<typename _MatrixType, int QRPreconditioner>
RealScalar Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_prescribedThreshold protected

template<typename _MatrixType, int QRPreconditioner>
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_qr_precond_morecols protected

template<typename _MatrixType, int QRPreconditioner>
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_qr_precond_morerows protected

template<typename _MatrixType, int QRPreconditioner>
Index Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_rows protected

template<typename _MatrixType, int QRPreconditioner>
MatrixType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_scaledMatrix protected

template<typename _MatrixType, int QRPreconditioner>
SingularValuesType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_singularValues protected

template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_usePrescribedThreshold protected

template<typename _MatrixType, int QRPreconditioner>
WorkMatrixType Eigen::JacobiSVD<_MatrixType, QRPreconditioner>::m_workMatrix protected