Eigen::JacobiSVD< MatrixType_, Options_ > Class Template Reference

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

+ Inheritance diagram for Eigen::JacobiSVD< MatrixType_, Options_ >:

Public Types

enum  {
  Options ,
  QRPreconditioner ,
  RowsAtCompileTime ,
  ColsAtCompileTime ,
  DiagSizeAtCompileTime ,
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime ,
  MaxDiagSizeAtCompileTime ,
  MatrixOptions
}
 
typedef Base::Index Index
 
typedef MatrixType_ MatrixType
 
typedef Base::MatrixUType MatrixUType
 
typedef Base::MatrixVType MatrixVType
 
typedef Base::RealScalar RealScalar
 
typedef Base::Scalar Scalar
 
typedef Base::SingularValuesType SingularValuesType
 
typedef Matrix< Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTimeWorkMatrixType
 
- Public Types inherited from Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >
enum  
 
typedef Eigen::Index Index
 
typedef internal::traits< JacobiSVD< MatrixType_, Options_ > >::MatrixType MatrixType
 
typedef internal::make_proper_matrix_type< Scalar, RowsAtCompileTime, MatrixUColsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MatrixUMaxColsAtCompileTime >::type MatrixUType
 
typedef internal::make_proper_matrix_type< Scalar, ColsAtCompileTime, MatrixVColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MatrixVMaxColsAtCompileTime >::type MatrixVType
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
 
typedef Eigen::internal::traits< SVDBase >::StorageIndex StorageIndex
 
- Public Types inherited from Eigen::SolverBase< Derived >
enum  {
  RowsAtCompileTime ,
  ColsAtCompileTime ,
  SizeAtCompileTime ,
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime ,
  MaxSizeAtCompileTime ,
  IsVectorAtCompileTime ,
  NumDimensions
}
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
typedef EigenBase< Derived > Base
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const Derived > ConstTransposeReturnType
 
typedef internal::traits< Derived >::Scalar Scalar
 
- Public Types inherited from Eigen::EigenBase< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 
typedef internal::traits< Derived >::StorageKind StorageKind
 

Public Member Functions

EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
JacobiSVDcompute (const MatrixType &matrix)
 Method performing the decomposition of given matrix. Computes Thin/Full unitaries U/V if specified using the Options template parameter or the class constructor. More...
 
EIGEN_DEPRECATED JacobiSVDcompute (const MatrixType &matrix, unsigned int computationOptions)
 Method performing the decomposition of given matrix, as specified by the computationOptions parameter. More...
 
 JacobiSVD ()
 Default Constructor. More...
 
 JacobiSVD (const MatrixType &matrix)
 Constructor performing the decomposition of given matrix, using the custom options specified with the Options template paramter. More...
 
 JacobiSVD (const MatrixType &matrix, unsigned int computationOptions)
 Constructor performing the decomposition of given matrix using specified options for computing unitaries. More...
 
 JacobiSVD (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
EIGEN_DEPRECATED JacobiSVD (Index rows, Index cols, unsigned int computationOptions)
 Default Constructor with memory preallocation. More...
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
- Public Member Functions inherited from Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >
Index cols () const
 
bool computeU () const
 
bool computeV () const
 
JacobiSVD< MatrixType_, Options_ > & derived ()
 
const JacobiSVD< MatrixType_, Options_ > & derived () const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
const MatrixUTypematrixU () const
 
const MatrixVTypematrixV () const
 
Index nonzeroSingularValues () const
 
Index rank () const
 
Index rows () const
 
JacobiSVD< MatrixType_, Options_ > & setThreshold (const RealScalar &threshold)
 
JacobiSVD< MatrixType_, Options_ > & setThreshold (Default_t)
 
const SingularValuesTypesingularValues () const
 
const Solve< JacobiSVD< MatrixType_, Options_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
- Public Member Functions inherited from Eigen::SolverBase< Derived >
const AdjointReturnType adjoint () const
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
 ~SolverBase ()
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
template<typename Dest >
void addTo (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheLeft (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheRight (Dest &dst) const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & const_cast_derived () const
 
const Derived & const_derived () const
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Dest >
void evalTo (Dest &dst) const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
template<typename Dest >
void subTo (Dest &dst) const
 

Protected Member Functions

 EIGEN_STATIC_ASSERT (!(ShouldComputeThinU &&int(QRPreconditioner)==int(FullPivHouseholderQRPreconditioner)) &&!(ShouldComputeThinU &&int(QRPreconditioner)==int(FullPivHouseholderQRPreconditioner)), "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. " "Use the ColPivHouseholderQR preconditioner instead.") template< typename MatrixType__
 
- Protected Member Functions inherited from Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >
void _check_compute_assertions () const
 
void _check_solve_assertion (const Rhs &b) const
 
bool allocate (Index rows, Index cols, unsigned int computationOptions)
 
 SVDBase ()
 Default Constructor. More...
 
- Protected Member Functions inherited from Eigen::SolverBase< Derived >
template<bool Transpose_, typename Rhs >
void _check_solve_assertion (const Rhs &b) const
 

Protected Attributes

internal::qr_preconditioner_impl< MatrixType, Options, QRPreconditioner, internal::PreconditionIfMoreColsThanRowsm_qr_precond_morecols
 
internal::qr_preconditioner_impl< MatrixType, Options, QRPreconditioner, internal::PreconditionIfMoreRowsThanColsm_qr_precond_morerows
 
MatrixType m_scaledMatrix
 
WorkMatrixType m_workMatrix
 
int Options__
 
- Protected Attributes inherited from Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >
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
 
Index m_rows
 
SingularValuesType m_singularValues
 
bool m_usePrescribedThreshold
 

Private Types

typedef SVDBase< JacobiSVDBase
 

Private Member Functions

void allocate (Index rows, Index cols, unsigned int computationOptions)
 
JacobiSVDcompute_impl (const MatrixType &matrix, unsigned int computationOptions)
 

Additional Inherited Members

- Static Public Attributes inherited from Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >
static constexpr bool ShouldComputeFullU
 
static constexpr bool ShouldComputeFullV
 
static constexpr bool ShouldComputeThinU
 
static constexpr bool ShouldComputeThinV
 

Detailed Description

template<typename MatrixType_, int Options_>
class Eigen::JacobiSVD< MatrixType_, Options_ >

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
Optionsthis optional parameter allows one to specify the type of QR decomposition that will be used internally for the R-SVD step for non-square matrices. Additionally, it allows one to specify whether to compute thin or full unitaries U and V. 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, ComputeThinU | ComputeThinV> svd(m);
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;
Matrix3f m
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< float, 3, 1 > Vector3f
3×1 vector of type float.
Definition: Matrix.h:501
Matrix< float, Dynamic, Dynamic > MatrixXf
Dynamic×Dynamic matrix of type float.
Definition: Matrix.h:501

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 QR preconditioners that can be set with Options template parameter 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.

One may also use the Options template parameter to specify how the unitaries should be computed. The options are ComputeThinU, ComputeThinV, ComputeFullU, ComputeFullV. It is not possible to request both the thin and full versions of a unitary. By default, unitaries will not be computed.

You can set the QRPreconditioner and unitary options together: JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner | ComputeThinU | ComputeFullV>

See also
MatrixBase::jacobiSvd()

Definition at line 514 of file JacobiSVD.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , int Options_>
typedef SVDBase<JacobiSVD> Eigen::JacobiSVD< MatrixType_, Options_ >::Base
private

Definition at line 515 of file JacobiSVD.h.

◆ Index

template<typename MatrixType_ , int Options_>
typedef Base::Index Eigen::JacobiSVD< MatrixType_, Options_ >::Index

Definition at line 521 of file JacobiSVD.h.

◆ MatrixType

template<typename MatrixType_ , int Options_>
typedef MatrixType_ Eigen::JacobiSVD< MatrixType_, Options_ >::MatrixType

Definition at line 518 of file JacobiSVD.h.

◆ MatrixUType

template<typename MatrixType_ , int Options_>
typedef Base::MatrixUType Eigen::JacobiSVD< MatrixType_, Options_ >::MatrixUType

Definition at line 534 of file JacobiSVD.h.

◆ MatrixVType

template<typename MatrixType_ , int Options_>
typedef Base::MatrixVType Eigen::JacobiSVD< MatrixType_, Options_ >::MatrixVType

Definition at line 535 of file JacobiSVD.h.

◆ RealScalar

template<typename MatrixType_ , int Options_>
typedef Base::RealScalar Eigen::JacobiSVD< MatrixType_, Options_ >::RealScalar

Definition at line 520 of file JacobiSVD.h.

◆ Scalar

template<typename MatrixType_ , int Options_>
typedef Base::Scalar Eigen::JacobiSVD< MatrixType_, Options_ >::Scalar

Definition at line 519 of file JacobiSVD.h.

◆ SingularValuesType

template<typename MatrixType_ , int Options_>
typedef Base::SingularValuesType Eigen::JacobiSVD< MatrixType_, Options_ >::SingularValuesType

Definition at line 536 of file JacobiSVD.h.

◆ WorkMatrixType

Definition at line 539 of file JacobiSVD.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , int Options_>
anonymous enum
Enumerator
Options 
QRPreconditioner 
RowsAtCompileTime 
ColsAtCompileTime 
DiagSizeAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
MaxDiagSizeAtCompileTime 
MatrixOptions 

Definition at line 522 of file JacobiSVD.h.

522  {
523  Options = Options_,
532  };
@ MaxDiagSizeAtCompileTime
Definition: JacobiSVD.h:530
@ ColsAtCompileTime
Definition: SVDBase.h:139
@ DiagSizeAtCompileTime
Definition: SVDBase.h:140
@ MaxRowsAtCompileTime
Definition: SVDBase.h:141
@ MaxDiagSizeAtCompileTime
Definition: SVDBase.h:143
@ MaxColsAtCompileTime
Definition: SVDBase.h:142
@ RowsAtCompileTime
Definition: SVDBase.h:138
constexpr int get_qr_preconditioner(int options)
Definition: SVDBase.h:31

Constructor & Destructor Documentation

◆ JacobiSVD() [1/5]

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::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 546 of file JacobiSVD.h.

546 {}

◆ JacobiSVD() [2/5]

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD ( Index  rows,
Index  cols 
)
inline

Default Constructor with memory preallocation.

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

See also
JacobiSVD()

Definition at line 555 of file JacobiSVD.h.

EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: EigenBase.h:65
void allocate(Index rows, Index cols, unsigned int computationOptions)
Definition: JacobiSVD.h:674
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: EigenBase.h:62
constexpr int get_computation_options(int options)
Definition: SVDBase.h:33

◆ JacobiSVD() [3/5]

template<typename MatrixType_ , int Options_>
EIGEN_DEPRECATED Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD ( Index  rows,
Index  cols,
unsigned int  computationOptions 
)
inline

Default Constructor with memory preallocation.

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

One cannot request unitaries using both the Options template parameter and the constructor. If possible, prefer using the Options template parameter.

Parameters
computationOptionsspecify whether to compute Thin/Full unitaries U/V
See also
JacobiSVD()
Deprecated:
Will be removed in the next major Eigen version. Options should be specified in the Options template parameter.

Definition at line 572 of file JacobiSVD.h.

572  {
573  internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, rows, cols);
574  allocate(rows, cols, computationOptions);
575  }

◆ JacobiSVD() [4/5]

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD ( const MatrixType matrix)
inlineexplicit

Constructor performing the decomposition of given matrix, using the custom options specified with the Options template paramter.

Parameters
matrixthe matrix to decompose

Definition at line 582 of file JacobiSVD.h.

JacobiSVD & compute_impl(const MatrixType &matrix, unsigned int computationOptions)
Definition: JacobiSVD.h:689

◆ JacobiSVD() [5/5]

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD ( const MatrixType matrix,
unsigned int  computationOptions 
)
inline

Constructor performing the decomposition of given matrix using specified options for computing unitaries.

One cannot request unitiaries using both the Options template parameter and the constructor. If possible, prefer using the Options template parameter.

Parameters
matrixthe matrix to decompose
computationOptionsspecify whether to compute Thin/Full unitaries U/V
Deprecated:
Will be removed in the next major Eigen version. Options should be specified in the Options template parameter.

Definition at line 597 of file JacobiSVD.h.

597  {
598  internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, matrix.rows(), matrix.cols());
599  compute_impl(matrix, computationOptions);
600  }

Member Function Documentation

◆ allocate()

template<typename MatrixType , int Options>
void Eigen::JacobiSVD< MatrixType, Options >::allocate ( Index  rows,
Index  cols,
unsigned int  computationOptions 
)
private

Definition at line 674 of file JacobiSVD.h.

674  {
675  if (Base::allocate(rows, cols, computationOptions)) return;
676 
679  "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
680  "Use the ColPivHouseholderQR preconditioner instead.");
681 
683  if(m_cols>m_rows) m_qr_precond_morecols.allocate(*this);
684  if(m_rows>m_cols) m_qr_precond_morerows.allocate(*this);
685  if(m_rows!=m_cols) m_scaledMatrix.resize(rows,cols);
686 }
#define eigen_assert(x)
Definition: Macros.h:902
internal::qr_preconditioner_impl< MatrixType, Options, QRPreconditioner, internal::PreconditionIfMoreColsThanRows > m_qr_precond_morecols
Definition: JacobiSVD.h:666
internal::qr_preconditioner_impl< MatrixType, Options, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols > m_qr_precond_morerows
Definition: JacobiSVD.h:668
MatrixType m_scaledMatrix
Definition: JacobiSVD.h:670
WorkMatrixType m_workMatrix
Definition: JacobiSVD.h:669
constexpr void resize(Index rows, Index cols)
bool allocate(Index rows, Index cols, unsigned int computationOptions)
Definition: SVDBase.h:410
@ FullPivHouseholderQRPreconditioner
Definition: Constants.h:433

◆ cols()

template<typename MatrixType_ , int Options_>
EIGEN_CONSTEXPR Index Eigen::EigenBase< Derived >::cols ( void  )
inline
Returns
the number of columns.
See also
rows(), ColsAtCompileTime

Definition at line 65 of file EigenBase.h.

65 { return derived().cols(); }
JacobiSVD< MatrixType_, Options_ > & derived()
Definition: SVDBase.h:163

◆ compute() [1/2]

template<typename MatrixType_ , int Options_>
JacobiSVD& Eigen::JacobiSVD< MatrixType_, Options_ >::compute ( const MatrixType matrix)
inline

Method performing the decomposition of given matrix. Computes Thin/Full unitaries U/V if specified using the Options template parameter or the class constructor.

Parameters
matrixthe matrix to decompose

Definition at line 607 of file JacobiSVD.h.

◆ compute() [2/2]

template<typename MatrixType_ , int Options_>
EIGEN_DEPRECATED JacobiSVD& Eigen::JacobiSVD< MatrixType_, Options_ >::compute ( const MatrixType matrix,
unsigned int  computationOptions 
)
inline

Method performing the decomposition of given matrix, as specified by the computationOptions parameter.

Parameters
matrixthe matrix to decompose
computationOptionsspecify whether to compute Thin/Full unitaries U/V
Deprecated:
Will be removed in the next major Eigen version. Options should be specified in the Options template parameter.

Definition at line 619 of file JacobiSVD.h.

619  {
620  internal::check_svd_options_assertions<MatrixType, Options>(m_computationOptions, matrix.rows(), matrix.cols());
621  return compute_impl(matrix, computationOptions);
622  }

◆ compute_impl()

template<typename MatrixType , int Options>
JacobiSVD< MatrixType, Options > & Eigen::JacobiSVD< MatrixType, Options >::compute_impl ( const MatrixType matrix,
unsigned int  computationOptions 
)
private

Definition at line 689 of file JacobiSVD.h.

690  {
691  using std::abs;
692 
693  allocate(matrix.rows(), matrix.cols(), computationOptions);
694 
695  // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
696  // only worsening the precision of U and V as we accumulate more rotations
697  const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
698 
699  // limit for denormal numbers to be considered zero in order to avoid infinite loops (see bug 286)
700  const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
701 
702  // Scaling factor to reduce over/under-flows
703  RealScalar scale = matrix.cwiseAbs().template maxCoeff<PropagateNaN>();
704  if (!(numext::isfinite)(scale)) {
705  m_isInitialized = true;
708  return *this;
709  }
710  if(numext::is_exactly_zero(scale)) scale = RealScalar(1);
711 
712 
714  if(m_rows!=m_cols)
715  {
716  m_scaledMatrix = matrix / scale;
719  }
720  else
721  {
722  m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize) / scale;
727  }
728 
729 
730  RealScalar maxDiagEntry = m_workMatrix.cwiseAbs().diagonal().maxCoeff();
731 
732  bool finished = false;
733  while(!finished)
734  {
735  finished = true;
736 
737  // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
738 
739  for(Index p = 1; p < m_diagSize; ++p)
740  {
741  for(Index q = 0; q < p; ++q)
742  {
743  // if this 2x2 sub-matrix is not diagonal already...
744  // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
745  // keep us iterating forever. Similarly, small denormal numbers are considered zero.
746  RealScalar threshold = numext::maxi<RealScalar>(considerAsZero, precision * maxDiagEntry);
748  {
749  finished = false;
750  // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
751  // the complex to real operation returns true if the updated 2x2 block is not already diagonal
752  if (internal::svd_precondition_2x2_block_to_be_real<MatrixType, Options>::run(m_workMatrix, *this, p, q,
753  maxDiagEntry)) {
754  JacobiRotation<RealScalar> j_left, j_right;
755  internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
756 
757  // accumulate resulting Jacobi rotations
758  m_workMatrix.applyOnTheLeft(p,q,j_left);
759  if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
760 
761  m_workMatrix.applyOnTheRight(p,q,j_right);
762  if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
763 
764  // keep track of the largest diagonal coefficient
765  maxDiagEntry = numext::maxi<RealScalar>(maxDiagEntry,numext::maxi<RealScalar>(abs(m_workMatrix.coeff(p,p)), abs(m_workMatrix.coeff(q,q))));
766  }
767  }
768  }
769  }
770  }
771 
772 
774  for(Index i = 0; i < m_diagSize; ++i)
775  {
776  // For a complex matrix, some diagonal coefficients might note have been
777  // treated by svd_precondition_2x2_block_to_be_real, and the imaginary part
778  // of some diagonal entry might not be null.
780  {
782  m_singularValues.coeffRef(i) = abs(a);
783  if(computeU()) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
784  }
785  else
786  {
787  // m_workMatrix.coeff(i,i) is already real, no difficulty:
789  m_singularValues.coeffRef(i) = abs(a);
790  if(computeU() && (a<RealScalar(0))) m_matrixU.col(i) = -m_matrixU.col(i);
791  }
792  }
793 
794  m_singularValues *= scale;
795 
796 
799  for(Index i = 0; i < m_diagSize; i++)
800  {
801  Index pos;
802  RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
803  if(numext::is_exactly_zero(maxRemainingSingularValue))
804  {
806  break;
807  }
808  if(pos)
809  {
810  pos += i;
811  std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
812  if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
813  if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
814  }
815  }
816 
817  m_isInitialized = true;
818  return *this;
819 }
const AbsReturnType abs() const
const ImagReturnType imag() const
RealReturnType real() const
float * p
void swap(const DenseBase< OtherDerived > &other)
Definition: DenseBase.h:418
Base::RealScalar RealScalar
Definition: JacobiSVD.h:520
Base::Index Index
Definition: JacobiSVD.h:521
Derived & setIdentity()
void applyOnTheRight(const EigenBase< OtherDerived > &other)
Definition: MatrixBase.h:532
constexpr const Scalar & coeff(Index rowId, Index colId) const
@ InvalidInput
Definition: Constants.h:453
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
void real_2x2_jacobi_svd(const MatrixType &matrix, Index p, Index q, JacobiRotation< RealScalar > *j_left, JacobiRotation< RealScalar > *j_right)
Definition: RealSvd2x2.h:21
void swap(scoped_array< T > &a, scoped_array< T > &b)
Definition: Memory.h:788
EIGEN_ALWAYS_INLINE bool() isfinite(const Eigen::bfloat16 &h)
Definition: BFloat16.h:790
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)

◆ EIGEN_STATIC_ASSERT()

template<typename MatrixType_ , int Options_>
Eigen::JacobiSVD< MatrixType_, Options_ >::EIGEN_STATIC_ASSERT ( ShouldComputeThinU &&int(QRPreconditioner)==int(FullPivHouseholderQRPreconditioner)) &&!(ShouldComputeThinU &&int(QRPreconditioner)==int(FullPivHouseholderQRPreconditioner),
"JacobiSVD< MatrixType_, Options_ >: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. " "Use the ColPivHouseholderQR preconditioner instead."   
)
protected

◆ rows()

template<typename MatrixType_ , int Options_>
EIGEN_CONSTEXPR Index Eigen::EigenBase< Derived >::rows ( void  )
inline
Returns
the number of rows.
See also
cols(), RowsAtCompileTime

Definition at line 62 of file EigenBase.h.

62 { return derived().rows(); }

Member Data Documentation

◆ m_qr_precond_morecols

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

Definition at line 666 of file JacobiSVD.h.

◆ m_qr_precond_morerows

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

Definition at line 668 of file JacobiSVD.h.

◆ m_scaledMatrix

template<typename MatrixType_ , int Options_>
MatrixType Eigen::JacobiSVD< MatrixType_, Options_ >::m_scaledMatrix
protected

Definition at line 670 of file JacobiSVD.h.

◆ m_workMatrix

template<typename MatrixType_ , int Options_>
WorkMatrixType Eigen::JacobiSVD< MatrixType_, Options_ >::m_workMatrix
protected

Definition at line 669 of file JacobiSVD.h.

◆ Options__

template<typename MatrixType_ , int Options_>
int Eigen::JacobiSVD< MatrixType_, Options_ >::Options__
protected

Definition at line 660 of file JacobiSVD.h.


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