Shapeworks Studio  2.1
Shape analysis software suite
List of all members | Public Types | Public Member Functions | Protected Attributes | Friends
Eigen::JacobiSVD< _MatrixType, QRPreconditioner > Class Template Reference

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

#include <JacobiSVD.h>

+ Collaboration diagram for Eigen::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 MatrixType::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
 
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 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...
 
const MatrixUTypematrixU () const
 
const MatrixVTypematrixV () const
 
const SingularValuesTypesingularValues () const
 
bool computeU () const
 
bool computeV () const
 
template<typename Rhs >
const internal::solve_retval< JacobiSVD, Rhs > solve (const MatrixBase< Rhs > &b) const
 
Index nonzeroSingularValues () const
 
Index rows () const
 
Index cols () const
 

Protected Attributes

MatrixUType m_matrixU
 
MatrixVType m_matrixV
 
SingularValuesType m_singularValues
 
WorkMatrixType m_workMatrix
 
bool m_isInitialized
 
bool m_isAllocated
 
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
 
internal::qr_preconditioner_impl< MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows > m_qr_precond_morecols
 
internal::qr_preconditioner_impl< MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols > m_qr_precond_morerows
 

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
 

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:

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 224 of file ForwardDeclarations.h.

Constructor & Destructor Documentation

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&).

Definition at line 528 of file JacobiSVD.h.

529  : m_isInitialized(false),
530  m_isAllocated(false),
531  m_computationOptions(0),
532  m_rows(-1), m_cols(-1)
533  {}
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()

Definition at line 542 of file JacobiSVD.h.

543  : m_isInitialized(false),
544  m_isAllocated(false),
545  m_computationOptions(0),
546  m_rows(-1), m_cols(-1)
547  {
548  allocate(rows, cols, computationOptions);
549  }
template<typename _MatrixType, int QRPreconditioner>
Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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.

Definition at line 561 of file JacobiSVD.h.

562  : m_isInitialized(false),
563  m_isAllocated(false),
564  m_computationOptions(0),
565  m_rows(-1), m_cols(-1)
566  {
567  compute(matrix, computationOptions);
568  }
Definition: math3d.h:219
JacobiSVD & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
Definition: JacobiSVD.h:742

Member Function Documentation

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.

Definition at line 742 of file JacobiSVD.h.

743 {
744  using std::abs;
745  allocate(matrix.rows(), matrix.cols(), computationOptions);
746 
747  // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
748  // only worsening the precision of U and V as we accumulate more rotations
749  const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
750 
751  // limit for very small denormal numbers to be considered zero in order to avoid infinite loops (see bug 286)
752  const RealScalar considerAsZero = RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min();
753 
754  /*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */
755 
756  if(!m_qr_precond_morecols.run(*this, matrix) && !m_qr_precond_morerows.run(*this, matrix))
757  {
758  m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize);
759  if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
760  if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
761  if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
762  if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
763  }
764 
765  /*** step 2. The main Jacobi SVD iteration. ***/
766 
767  bool finished = false;
768  while(!finished)
769  {
770  finished = true;
771 
772  // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
773 
774  for(Index p = 1; p < m_diagSize; ++p)
775  {
776  for(Index q = 0; q < p; ++q)
777  {
778  // if this 2x2 sub-matrix is not diagonal already...
779  // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
780  // keep us iterating forever. Similarly, small denormal numbers are considered zero.
781  using std::max;
782  RealScalar threshold = (max)(considerAsZero, precision * (max)(abs(m_workMatrix.coeff(p,p)),
783  abs(m_workMatrix.coeff(q,q))));
784  if((max)(abs(m_workMatrix.coeff(p,q)),abs(m_workMatrix.coeff(q,p))) > threshold)
785  {
786  finished = false;
787 
788  // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
789  internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q);
790  JacobiRotation<RealScalar> j_left, j_right;
791  internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
792 
793  // accumulate resulting Jacobi rotations
794  m_workMatrix.applyOnTheLeft(p,q,j_left);
795  if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
796 
797  m_workMatrix.applyOnTheRight(p,q,j_right);
798  if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
799  }
800  }
801  }
802  }
803 
804  /*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/
805 
806  for(Index i = 0; i < m_diagSize; ++i)
807  {
808  RealScalar a = abs(m_workMatrix.coeff(i,i));
809  m_singularValues.coeffRef(i) = a;
810  if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
811  }
812 
813  /*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/
814 
815  m_nonzeroSingularValues = m_diagSize;
816  for(Index i = 0; i < m_diagSize; i++)
817  {
818  Index pos;
819  RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
820  if(maxRemainingSingularValue == RealScalar(0))
821  {
822  m_nonzeroSingularValues = i;
823  break;
824  }
825  if(pos)
826  {
827  pos += i;
828  std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
829  if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
830  if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
831  }
832  }
833 
834  m_isInitialized = true;
835  return *this;
836 }
bool computeV() const
Definition: JacobiSVD.h:639
Definition: math3d.h:219
bool computeU() const
Definition: JacobiSVD.h:637
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).

Definition at line 588 of file JacobiSVD.h.

589  {
590  return compute(matrix, m_computationOptions);
591  }
Definition: math3d.h:219
JacobiSVD & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
Definition: JacobiSVD.h:742
template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::computeU ( ) const
inline
Returns
true if U (full or thin) is asked for in this SVD decomposition

Definition at line 637 of file JacobiSVD.h.

637 { return m_computeFullU || m_computeThinU; }
template<typename _MatrixType, int QRPreconditioner>
bool Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::computeV ( ) const
inline
Returns
true if V (full or thin) is asked for in this SVD decomposition

Definition at line 639 of file JacobiSVD.h.

639 { return m_computeFullV || m_computeThinV; }
template<typename _MatrixType, int QRPreconditioner>
const MatrixUType& Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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.

Definition at line 602 of file JacobiSVD.h.

603  {
604  eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
605  eigen_assert(computeU() && "This JacobiSVD decomposition didn't compute U. Did you ask for it?");
606  return m_matrixU;
607  }
bool computeU() const
Definition: JacobiSVD.h:637
template<typename _MatrixType, int QRPreconditioner>
const MatrixVType& Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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.

Definition at line 618 of file JacobiSVD.h.

619  {
620  eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
621  eigen_assert(computeV() && "This JacobiSVD decomposition didn't compute V. Did you ask for it?");
622  return m_matrixV;
623  }
bool computeV() const
Definition: JacobiSVD.h:639
template<typename _MatrixType, int QRPreconditioner>
Index Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::nonzeroSingularValues ( ) const
inline
Returns
the number of singular values that are not exactly 0

Definition at line 660 of file JacobiSVD.h.

661  {
662  eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
663  return m_nonzeroSingularValues;
664  }
template<typename _MatrixType, int QRPreconditioner>
const SingularValuesType& Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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.

Definition at line 630 of file JacobiSVD.h.

631  {
632  eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
633  return m_singularValues;
634  }
template<typename _MatrixType, int QRPreconditioner>
template<typename Rhs >
const internal::solve_retval<JacobiSVD, Rhs> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::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 $.

Definition at line 652 of file JacobiSVD.h.

653  {
654  eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
655  eigen_assert(computeU() && computeV() && "JacobiSVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
656  return internal::solve_retval<JacobiSVD, Rhs>(*this, b.derived());
657  }
bool computeV() const
Definition: JacobiSVD.h:639
bool computeU() const
Definition: JacobiSVD.h:637

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