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

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

#include <JacobiSVD.h>

Inheritance diagram for Eigen::JacobiSVD< _MatrixType, QRPreconditioner >:
Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >

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 _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef Base::MatrixUType MatrixUType
 
typedef Base::MatrixVType MatrixVType
 
typedef Base::SingularValuesType SingularValuesType
 
typedef internal::plain_row_type< MatrixType >::type RowType
 
typedef internal::plain_col_type< MatrixType >::type ColType
 
typedef Matrix< Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime > WorkMatrixType
 
- Public Types inherited from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >
enum  
 
typedef internal::traits< JacobiSVD< _MatrixType, QRPreconditioner > >::MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Eigen::Index Index
 
typedef Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime > MatrixUType
 
typedef Matrix< Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime > MatrixVType
 
typedef internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
 

Public Member Functions

 JacobiSVD ()
 Default Constructor. More...
 
 JacobiSVD (Index rows, Index cols, unsigned int computationOptions=0)
 Default Constructor with memory preallocation. More...
 
 JacobiSVD (const MatrixType &matrix, unsigned int computationOptions=0)
 Constructor performing the decomposition of given matrix. More...
 
JacobiSVDcompute (const MatrixType &matrix, unsigned int computationOptions)
 Method performing the decomposition of given matrix using custom options. More...
 
JacobiSVDcompute (const MatrixType &matrix)
 Method performing the decomposition of given matrix using current options. More...
 
- Public Member Functions inherited from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >
JacobiSVD< _MatrixType, QRPreconditioner > & derived ()
 
const JacobiSVD< _MatrixType, QRPreconditioner > & derived () const
 
const MatrixUTypematrixU () const
 
const MatrixVTypematrixV () const
 
const SingularValuesTypesingularValues () const
 
Index nonzeroSingularValues () const
 
Index rank () const
 
JacobiSVD< _MatrixType, QRPreconditioner > & setThreshold (const RealScalar &threshold)
 Allows to prescribe a threshold to be used by certain methods, such as rank() and solve(), which need to determine when singular values are to be considered nonzero. More...
 
JacobiSVD< _MatrixType, QRPreconditioner > & setThreshold (Default_t)
 Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold. More...
 
RealScalar threshold () const
 Returns the threshold that will be used by certain methods such as rank(). More...
 
bool computeU () const
 
bool computeV () const
 
Index rows () const
 
Index cols () const
 
const Solve< JacobiSVD< _MatrixType, QRPreconditioner >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
EIGEN_DEVICE_FUNC void _solve_impl (const RhsType &rhs, DstType &dst) const
 
void _solve_impl (const RhsType &rhs, DstType &dst) const
 

Protected Attributes

WorkMatrixType m_workMatrix
 
internal::qr_preconditioner_impl< MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows > m_qr_precond_morecols
 
internal::qr_preconditioner_impl< MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols > m_qr_precond_morerows
 
MatrixType m_scaledMatrix
 
- Protected Attributes inherited from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >
MatrixUType m_matrixU
 
MatrixVType m_matrixV
 
SingularValuesType m_singularValues
 
bool m_isInitialized
 
bool m_isAllocated
 
bool m_usePrescribedThreshold
 
bool m_computeFullU
 
bool m_computeThinU
 
bool m_computeFullV
 
bool m_computeThinV
 
unsigned int m_computationOptions
 
Index m_nonzeroSingularValues
 
Index m_rows
 
Index m_cols
 
Index m_diagSize
 
RealScalar m_prescribedThreshold
 

Friends

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

Additional Inherited Members

- Protected Member Functions inherited from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >
bool allocate (Index rows, Index cols, unsigned int computationOptions)
 
 SVDBase ()
 Default Constructor. More...
 
- Static Protected Member Functions inherited from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >
static void check_template_parameters ()
 

Detailed Description

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

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:

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

Constructor & Destructor Documentation

§ JacobiSVD() [1/3]

template<typename _MatrixType, int QRPreconditioner>
Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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/3]

template<typename _MatrixType, int QRPreconditioner>
Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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() [3/3]

template<typename _MatrixType, int QRPreconditioner>
Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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.

Member Function Documentation

§ compute() [1/2]

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

§ compute() [2/2]

template<typename _MatrixType, int QRPreconditioner>
JacobiSVD& Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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).


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