Public Types | Public Member Functions | Protected Attributes | Private Types | Private Member Functions | Friends | List of all members
Eigen::JacobiSVD Class Reference

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

#include <ForwardDeclarations.h>

Public Types

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

Public Member Functions

Index cols () const
 
JacobiSVDcompute (const MatrixType &matrix)
 Method performing the decomposition of given matrix using current options. More...
 
JacobiSVDcompute (const MatrixType &matrix, unsigned int computationOptions)
 Method performing the decomposition of given matrix using custom options. More...
 
bool computeU () const
 
bool computeV () const
 
 JacobiSVD ()
 Default Constructor. More...
 
 JacobiSVD (const MatrixType &matrix, unsigned int computationOptions=0)
 Constructor performing the decomposition of given matrix. More...
 
 JacobiSVD (Index rows, Index cols, unsigned int computationOptions=0)
 Default Constructor with memory preallocation. More...
 
Index rank () const
 
Index rows () 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
 
ComputationInfo m_info
 
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::PreconditionIfMoreColsThanRowsm_qr_precond_morecols
 
internal::qr_preconditioner_impl< MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanColsm_qr_precond_morerows
 
Index m_rows
 
MatrixType m_scaledMatrix
 
SingularValuesType m_singularValues
 
bool m_usePrescribedThreshold
 
WorkMatrixType m_workMatrix
 

Private Types

typedef SVDBase< JacobiSVDBase
 

Private Member Functions

void allocate (Index rows, Index cols, unsigned int computationOptions)
 

Friends

template<typename __MatrixType , int _QRPreconditioner, int _Case, bool _DoAnything>
struct internal::qr_preconditioner_impl
 
template<typename __MatrixType , int _QRPreconditioner, bool _IsComplex>
struct internal::svd_precondition_2x2_block_to_be_real
 

Detailed Description

Two-sided Jacobi SVD decomposition of a rectangular matrix.

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

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:

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

Definition at line 278 of file ForwardDeclarations.h.

Member Typedef Documentation

◆ Base

Definition at line 491 of file JacobiSVD.h.

◆ ColType

Definition at line 512 of file JacobiSVD.h.

◆ MatrixType

typedef _MatrixType Eigen::JacobiSVD::MatrixType

Definition at line 494 of file JacobiSVD.h.

◆ MatrixUType

Definition at line 507 of file JacobiSVD.h.

◆ MatrixVType

Definition at line 508 of file JacobiSVD.h.

◆ RealScalar

Definition at line 496 of file JacobiSVD.h.

◆ RowType

Definition at line 511 of file JacobiSVD.h.

◆ Scalar

typedef MatrixType::Scalar Eigen::JacobiSVD::Scalar

Definition at line 495 of file JacobiSVD.h.

◆ SingularValuesType

Definition at line 509 of file JacobiSVD.h.

◆ WorkMatrixType

Definition at line 515 of file JacobiSVD.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
DiagSizeAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
MaxDiagSizeAtCompileTime 
MatrixOptions 

Definition at line 497 of file JacobiSVD.h.

Constructor & Destructor Documentation

◆ JacobiSVD() [1/3]

Eigen::JacobiSVD::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 522 of file JacobiSVD.h.

◆ JacobiSVD() [2/3]

Eigen::JacobiSVD::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()

Definition at line 532 of file JacobiSVD.h.

◆ JacobiSVD() [3/3]

Eigen::JacobiSVD::JacobiSVD ( const MatrixType matrix,
unsigned int  computationOptions = 0 
)
inlineexplicit

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.

Definition at line 547 of file JacobiSVD.h.

Member Function Documentation

◆ allocate()

void Eigen::JacobiSVD::allocate ( Eigen::Index  rows,
Eigen::Index  cols,
unsigned int  computationOptions 
)
private

Definition at line 615 of file JacobiSVD.h.

◆ cols()

Index Eigen::SVDBase::cols
inline

Definition at line 213 of file SVDBase.h.

◆ compute() [1/2]

JacobiSVD& Eigen::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).

Definition at line 570 of file JacobiSVD.h.

◆ compute() [2/2]

JacobiSVD< MatrixType, QRPreconditioner > & Eigen::JacobiSVD::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.

Definition at line 666 of file JacobiSVD.h.

◆ computeU()

bool Eigen::SVDBase::computeU
inline
Returns
true if U (full or thin) is asked for in this SVD decomposition

Definition at line 208 of file SVDBase.h.

◆ computeV()

bool Eigen::SVDBase::computeV
inline
Returns
true if V (full or thin) is asked for in this SVD decomposition

Definition at line 210 of file SVDBase.h.

◆ rank()

Index Eigen::SVDBase::rank
inline
Returns
the rank of the matrix of which *this is the SVD.
Note
This method has to determine which singular values should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 148 of file SVDBase.h.

◆ rows()

Index Eigen::SVDBase::rows
inline

Definition at line 212 of file SVDBase.h.

Friends And Related Function Documentation

◆ internal::qr_preconditioner_impl

template<typename __MatrixType , int _QRPreconditioner, int _Case, bool _DoAnything>
friend struct internal::qr_preconditioner_impl
friend

Definition at line 607 of file JacobiSVD.h.

◆ internal::svd_precondition_2x2_block_to_be_real

template<typename __MatrixType , int _QRPreconditioner, bool _IsComplex>
friend struct internal::svd_precondition_2x2_block_to_be_real
friend

Definition at line 605 of file JacobiSVD.h.

Member Data Documentation

◆ m_cols

Index Eigen::SVDBase::m_cols
protected

Definition at line 280 of file SVDBase.h.

◆ m_computationOptions

unsigned int Eigen::SVDBase::m_computationOptions
protected

Definition at line 279 of file SVDBase.h.

◆ m_computeFullU

bool Eigen::SVDBase::m_computeFullU
protected

Definition at line 277 of file SVDBase.h.

◆ m_computeFullV

bool Eigen::SVDBase::m_computeFullV
protected

Definition at line 278 of file SVDBase.h.

◆ m_computeThinU

bool Eigen::SVDBase::m_computeThinU
protected

Definition at line 277 of file SVDBase.h.

◆ m_computeThinV

bool Eigen::SVDBase::m_computeThinV
protected

Definition at line 278 of file SVDBase.h.

◆ m_diagSize

Index Eigen::SVDBase::m_diagSize
protected

Definition at line 280 of file SVDBase.h.

◆ m_info

ComputationInfo Eigen::SVDBase::m_info
protected

Definition at line 275 of file SVDBase.h.

◆ m_isAllocated

bool Eigen::SVDBase::m_isAllocated
protected

Definition at line 276 of file SVDBase.h.

◆ m_isInitialized

bool Eigen::SVDBase::m_isInitialized
protected

Definition at line 276 of file SVDBase.h.

◆ m_matrixU

MatrixUType Eigen::SVDBase::m_matrixU
protected

Definition at line 272 of file SVDBase.h.

◆ m_matrixV

MatrixVType Eigen::SVDBase::m_matrixV
protected

Definition at line 273 of file SVDBase.h.

◆ m_nonzeroSingularValues

Index Eigen::SVDBase::m_nonzeroSingularValues
protected

Definition at line 280 of file SVDBase.h.

◆ m_prescribedThreshold

RealScalar Eigen::SVDBase::m_prescribedThreshold
protected

Definition at line 281 of file SVDBase.h.

◆ m_qr_precond_morecols

internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> Eigen::JacobiSVD::m_qr_precond_morecols
protected

Definition at line 609 of file JacobiSVD.h.

◆ m_qr_precond_morerows

internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> Eigen::JacobiSVD::m_qr_precond_morerows
protected

Definition at line 610 of file JacobiSVD.h.

◆ m_rows

Index Eigen::SVDBase::m_rows
protected

Definition at line 280 of file SVDBase.h.

◆ m_scaledMatrix

MatrixType Eigen::JacobiSVD::m_scaledMatrix
protected

Definition at line 611 of file JacobiSVD.h.

◆ m_singularValues

SingularValuesType Eigen::SVDBase::m_singularValues
protected

Definition at line 274 of file SVDBase.h.

◆ m_usePrescribedThreshold

bool Eigen::SVDBase::m_usePrescribedThreshold
protected

Definition at line 276 of file SVDBase.h.

◆ m_workMatrix

WorkMatrixType Eigen::JacobiSVD::m_workMatrix
protected

Definition at line 602 of file JacobiSVD.h.


The documentation for this class was generated from the following files:
svd
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf > svd(m, ComputeThinU|ComputeThinV)
Eigen::ComputeThinU
@ ComputeThinU
Definition: Constants.h:395
m
Matrix3f m
Definition: AngleAxis_mimic_euler.cpp:1
Eigen::ComputeThinV
@ ComputeThinV
Definition: Constants.h:399


gtsam
Author(s):
autogenerated on Sat Jun 1 2024 03:09:39