Public Types | Public Member Functions | Protected Attributes | List of all members
JacobiSVD< _MatrixType, QRPreconditioner > Class Template Reference

Two-sided Jacobi SVD decomposition of a rectangular matrix. More...

#include <JacobiSVD.h>

Public Types

enum  {
  RowsAtCompileTime,
  ColsAtCompileTime,
  DiagSizeAtCompileTime,
  MaxRowsAtCompileTime,
  MaxColsAtCompileTime,
  MaxDiagSizeAtCompileTime,
  MatrixOptions
}
typedef
internal::plain_col_type
< MatrixType >::type 
ColType
typedef MatrixType::Index Index
typedef _MatrixType MatrixType
typedef Matrix< Scalar,
RowsAtCompileTime,
RowsAtCompileTime,
MatrixOptions,
MaxRowsAtCompileTime,
MaxRowsAtCompileTime
MatrixUType
typedef Matrix< Scalar,
ColsAtCompileTime,
ColsAtCompileTime,
MatrixOptions,
MaxColsAtCompileTime,
MaxColsAtCompileTime
MatrixVType
typedef NumTraits< typename
MatrixType::Scalar >::Real 
RealScalar
typedef
internal::plain_row_type
< MatrixType >::type 
RowType
typedef MatrixType::Scalar Scalar
typedef
internal::plain_diag_type
< MatrixType, RealScalar >
::type 
SingularValuesType
typedef Matrix< Scalar,
DiagSizeAtCompileTime,
DiagSizeAtCompileTime,
MatrixOptions,
MaxDiagSizeAtCompileTime,
MaxDiagSizeAtCompileTime
WorkMatrixType

Public Member Functions

Index cols () const
JacobiSVDcompute (const MatrixType &matrix, unsigned int computationOptions)
 Method performing the decomposition of given matrix using custom options.
JacobiSVDcompute (const MatrixType &matrix)
 Method performing the decomposition of given matrix using current options.
bool computeU () const
bool computeV () const
 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.
const MatrixUTypematrixU () const
const MatrixVTypematrixV () const
Index nonzeroSingularValues () const
Index rows () const
const SingularValuesTypesingularValues () const
template<typename Rhs >
const internal::solve_retval
< JacobiSVD, Rhs > 
solve (const MatrixBase< Rhs > &b) const

Protected Attributes

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
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
SingularValuesType m_singularValues
WorkMatrixType m_workMatrix

Detailed Description

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

Two-sided Jacobi SVD decomposition of a rectangular matrix.

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

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:

See Also
MatrixBase::jacobiSvd()

Member Typedef Documentation

typedef internal::plain_col_type<MatrixType>::type ColType
typedef MatrixType::Index Index
typedef _MatrixType MatrixType
typedef NumTraits<typename MatrixType::Scalar>::Real RealScalar
typedef internal::plain_row_type<MatrixType>::type RowType
typedef MatrixType::Scalar Scalar
typedef internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType

Member Enumeration Documentation

anonymous enum
Enumerator:
RowsAtCompileTime 
ColsAtCompileTime 
DiagSizeAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
MaxDiagSizeAtCompileTime 
MatrixOptions 

Constructor & Destructor Documentation

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

See Also
JacobiSVD()
JacobiSVD ( const MatrixType matrix,
unsigned int  computationOptions = 0 
)
inline

Constructor performing the decomposition of given matrix.

Parameters
matrixthe matrix to decompose
computationOptionsoptional 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.

References JacobiSVD< _MatrixType, QRPreconditioner >::compute().

Member Function Documentation

Index cols ( void  ) const
inline
JacobiSVD< MatrixType, QRPreconditioner > & compute ( const MatrixType matrix,
unsigned int  computationOptions 
)

Method performing the decomposition of given matrix using custom options.

Parameters
matrixthe matrix to decompose
computationOptionsoptional 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.

References abs(), p, q, Eigen::internal::real_2x2_jacobi_svd(), and JacobiRotation< Scalar >::transpose().

Referenced by JacobiSVD< _MatrixType, QRPreconditioner >::compute(), and JacobiSVD< _MatrixType, QRPreconditioner >::JacobiSVD().

JacobiSVD& compute ( const MatrixType matrix)
inline

Method performing the decomposition of given matrix using current options.

Parameters
matrixthe matrix to decompose

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

References JacobiSVD< _MatrixType, QRPreconditioner >::compute(), and JacobiSVD< _MatrixType, QRPreconditioner >::m_computationOptions.

bool computeU ( ) const
inline
bool computeV ( ) const
inline
const MatrixUType& matrixU ( ) const
inline
Returns
the U matrix.

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.

References JacobiSVD< _MatrixType, QRPreconditioner >::computeU(), eigen_assert, JacobiSVD< _MatrixType, QRPreconditioner >::m_isInitialized, and JacobiSVD< _MatrixType, QRPreconditioner >::m_matrixU.

Referenced by Transform< _Scalar, _Dim, _Mode, _Options >::computeRotationScaling(), Transform< _Scalar, _Dim, _Mode, _Options >::computeScalingRotation(), and Eigen::umeyama().

const MatrixVType& matrixV ( ) const
inline
Returns
the V matrix.

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.

References JacobiSVD< _MatrixType, QRPreconditioner >::computeV(), eigen_assert, JacobiSVD< _MatrixType, QRPreconditioner >::m_isInitialized, and JacobiSVD< _MatrixType, QRPreconditioner >::m_matrixV.

Referenced by Transform< _Scalar, _Dim, _Mode, _Options >::computeRotationScaling(), Transform< _Scalar, _Dim, _Mode, _Options >::computeScalingRotation(), QuaternionBase< Derived >::setFromTwoVectors(), and Eigen::umeyama().

Index nonzeroSingularValues ( ) const
inline
Index rows ( void  ) const
inline
const SingularValuesType& singularValues ( ) const
inline
Returns
the vector of singular values.

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.

References eigen_assert, JacobiSVD< _MatrixType, QRPreconditioner >::m_isInitialized, and JacobiSVD< _MatrixType, QRPreconditioner >::m_singularValues.

Referenced by Transform< _Scalar, _Dim, _Mode, _Options >::computeRotationScaling(), Transform< _Scalar, _Dim, _Mode, _Options >::computeScalingRotation(), and Eigen::umeyama().

const internal::solve_retval<JacobiSVD, Rhs> solve ( const MatrixBase< Rhs > &  b) const
inline
Returns
a (least squares) solution of $ A x = b $ using the current SVD decomposition of A.
Parameters
bthe right-hand-side of the equation to solve.
Note
Solving requires both U and V to be computed. Thin U and V are enough, there is no need for full U or V.
SVD solving is implicitly least-squares. Thus, this method serves both purposes of exact solving and least-squares solving. In other words, the returned solution is guaranteed to minimize the Euclidean norm $ \Vert A x - b \Vert $.

References JacobiSVD< _MatrixType, QRPreconditioner >::computeU(), JacobiSVD< _MatrixType, QRPreconditioner >::computeV(), eigen_assert, and JacobiSVD< _MatrixType, QRPreconditioner >::m_isInitialized.

Member Data Documentation

Index m_cols
protected
unsigned int m_computationOptions
protected
bool m_computeFullU
protected
bool m_computeFullV
protected
bool m_computeThinU
protected
bool m_computeThinV
protected
Index m_diagSize
protected
bool m_isAllocated
protected
bool m_isInitialized
protected
MatrixUType m_matrixU
protected
MatrixVType m_matrixV
protected
Index m_nonzeroSingularValues
protected
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols
protected
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows
protected
Index m_rows
protected
SingularValuesType m_singularValues
protected
WorkMatrixType m_workMatrix
protected

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