Eigen::SVDBase< Derived > Class Template Reference

Base class of SVD algorithms. More...

#include <SVDBase.h>

+ Inheritance diagram for Eigen::SVDBase< Derived >:

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime , ColsAtCompileTime = MatrixType::ColsAtCompileTime , DiagSizeAtCompileTime = internal::min_size_prefer_dynamic(RowsAtCompileTime, ColsAtCompileTime) , MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime ,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime , MaxDiagSizeAtCompileTime = internal::min_size_prefer_fixed(MaxRowsAtCompileTime, MaxColsAtCompileTime) , MatrixOptions = internal::traits<MatrixType>::Options , MatrixUColsAtCompileTime = internal::traits<Derived>::MatrixUColsAtCompileTime ,
  MatrixVColsAtCompileTime = internal::traits<Derived>::MatrixVColsAtCompileTime , MatrixUMaxColsAtCompileTime = internal::traits<Derived>::MatrixUMaxColsAtCompileTime , MatrixVMaxColsAtCompileTime = internal::traits<Derived>::MatrixVMaxColsAtCompileTime
}
 
typedef internal::traits< Derived >::MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef Eigen::internal::traits< SVDBase >::StorageIndex StorageIndex
 
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 internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
 
- Public Types inherited from Eigen::SolverBase< SVDBase< Derived > >
enum  
 
typedef EigenBase< SVDBase< Derived > > Base
 
typedef internal::traits< SVDBase< Derived > >::Scalar Scalar
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const SVDBase< Derived > > ConstTransposeReturnType
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
- 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

Derived & derived ()
 
const Derived & derived () const
 
const MatrixUTypematrixU () const
 
const MatrixVTypematrixV () const
 
const SingularValuesTypesingularValues () const
 
Index nonzeroSingularValues () const
 
Index rank () const
 
Derived & setThreshold (const RealScalar &threshold)
 
Derived & setThreshold (Default_t)
 
RealScalar threshold () const
 
bool computeU () const
 
bool computeV () const
 
Index rows () const
 
Index cols () const
 
Index diagSize () const
 
EIGEN_DEVICE_FUNC ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<bool Conjugate, typename RhsType , typename DstType >
void _solve_impl_transposed (const RhsType &rhs, DstType &dst) const
 
- Public Member Functions inherited from Eigen::SolverBase< SVDBase< Derived > >
 SolverBase ()
 
 ~SolverBase ()
 
const Solve< SVDBase< Derived >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const ConstTransposeReturnType transpose () const
 
const AdjointReturnType adjoint () const
 
constexpr EIGEN_DEVICE_FUNC SVDBase< Derived > & derived ()
 
constexpr EIGEN_DEVICE_FUNC const SVDBase< Derived > & derived () const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
constexpr EIGEN_DEVICE_FUNC Derived & derived ()
 
constexpr EIGEN_DEVICE_FUNC const Derived & derived () const
 
EIGEN_DEVICE_FUNC Derived & const_cast_derived () const
 
EIGEN_DEVICE_FUNC const Derived & const_derived () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
template<typename Dest >
EIGEN_DEVICE_FUNC void evalTo (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void addTo (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void subTo (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void applyThisOnTheRight (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void applyThisOnTheLeft (Dest &dst) const
 
template<typename Device >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DeviceWrapper< Derived, Device > device (Device &device)
 
template<typename Device >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DeviceWrapper< const Derived, Device > device (Device &device) const
 

Static Public Attributes

static constexpr bool ShouldComputeFullU = internal::traits<Derived>::ShouldComputeFullU
 
static constexpr bool ShouldComputeThinU = internal::traits<Derived>::ShouldComputeThinU
 
static constexpr bool ShouldComputeFullV = internal::traits<Derived>::ShouldComputeFullV
 
static constexpr bool ShouldComputeThinV = internal::traits<Derived>::ShouldComputeThinV
 

Protected Member Functions

void _check_compute_assertions () const
 
template<bool Transpose_, typename Rhs >
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< SVDBase< Derived > >
void _check_solve_assertion (const Rhs &b) const
 

Protected Attributes

MatrixUType m_matrixU
 
MatrixVType m_matrixV
 
SingularValuesType m_singularValues
 
ComputationInfo m_info
 
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
 
internal::variable_if_dynamic< Index, RowsAtCompileTimem_rows
 
internal::variable_if_dynamic< Index, ColsAtCompileTimem_cols
 
internal::variable_if_dynamic< Index, DiagSizeAtCompileTimem_diagSize
 
RealScalar m_prescribedThreshold
 

Friends

template<typename Derived_ >
struct internal::solve_assertion
 

Detailed Description

template<typename Derived>
class Eigen::SVDBase< Derived >

Base class of SVD algorithms.

Template Parameters
Derivedthe type of the actual SVD decomposition

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.

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.

The status of the computation can be retrieved using the info() method. Unless info() returns Success, the results should be not considered well defined.

If the input matrix has inf or nan coefficients, the result of the computation is undefined, and info() will return InvalidInput, but the computation is guaranteed to terminate in finite (and reasonable) time.

See also
class BDCSVD, class JacobiSVD

Member Typedef Documentation

◆ MatrixType

template<typename Derived >
typedef internal::traits<Derived>::MatrixType Eigen::SVDBase< Derived >::MatrixType

◆ MatrixUType

◆ MatrixVType

◆ RealScalar

template<typename Derived >
typedef NumTraits<typename MatrixType::Scalar>::Real Eigen::SVDBase< Derived >::RealScalar

◆ Scalar

template<typename Derived >
typedef MatrixType::Scalar Eigen::SVDBase< Derived >::Scalar

◆ SingularValuesType

template<typename Derived >
typedef internal::plain_diag_type<MatrixType, RealScalar>::type Eigen::SVDBase< Derived >::SingularValuesType

◆ StorageIndex

template<typename Derived >
typedef Eigen::internal::traits<SVDBase>::StorageIndex Eigen::SVDBase< Derived >::StorageIndex

Member Enumeration Documentation

◆ anonymous enum

template<typename Derived >
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
DiagSizeAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
MaxDiagSizeAtCompileTime 
MatrixOptions 
MatrixUColsAtCompileTime 
MatrixVColsAtCompileTime 
MatrixUMaxColsAtCompileTime 
MatrixVMaxColsAtCompileTime 
134  {
135  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
136  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
138  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
139  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
141  MatrixOptions = internal::traits<MatrixType>::Options,
142  MatrixUColsAtCompileTime = internal::traits<Derived>::MatrixUColsAtCompileTime,
143  MatrixVColsAtCompileTime = internal::traits<Derived>::MatrixVColsAtCompileTime,
144  MatrixUMaxColsAtCompileTime = internal::traits<Derived>::MatrixUMaxColsAtCompileTime,
145  MatrixVMaxColsAtCompileTime = internal::traits<Derived>::MatrixVMaxColsAtCompileTime
146  };
@ MatrixVColsAtCompileTime
Definition: SVDBase.h:143
@ MatrixVMaxColsAtCompileTime
Definition: SVDBase.h:145
@ MatrixOptions
Definition: SVDBase.h:141
@ ColsAtCompileTime
Definition: SVDBase.h:136
@ DiagSizeAtCompileTime
Definition: SVDBase.h:137
@ MaxRowsAtCompileTime
Definition: SVDBase.h:138
@ MaxDiagSizeAtCompileTime
Definition: SVDBase.h:140
@ MatrixUColsAtCompileTime
Definition: SVDBase.h:142
@ MaxColsAtCompileTime
Definition: SVDBase.h:139
@ RowsAtCompileTime
Definition: SVDBase.h:135
@ MatrixUMaxColsAtCompileTime
Definition: SVDBase.h:144
constexpr int min_size_prefer_fixed(A a, B b)
Definition: Meta.h:683
constexpr int min_size_prefer_dynamic(A a, B b)
Definition: Meta.h:668

Constructor & Destructor Documentation

◆ SVDBase()

template<typename Derived >
Eigen::SVDBase< Derived >::SVDBase ( )
inlineprotected

Default Constructor.

Default constructor of SVDBase

350  : m_matrixU(MatrixUType()),
353  m_info(Success),
354  m_isInitialized(false),
355  m_isAllocated(false),
361  m_computationOptions(internal::traits<Derived>::Options),
static constexpr bool ShouldComputeFullV
Definition: SVDBase.h:131
internal::make_proper_matrix_type< Scalar, RowsAtCompileTime, MatrixUColsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MatrixUMaxColsAtCompileTime >::type MatrixUType
Definition: SVDBase.h:153
bool m_usePrescribedThreshold
Definition: SVDBase.h:335
bool m_computeFullV
Definition: SVDBase.h:337
internal::variable_if_dynamic< Index, DiagSizeAtCompileTime > m_diagSize
Definition: SVDBase.h:342
static constexpr bool ShouldComputeThinV
Definition: SVDBase.h:132
ComputationInfo m_info
Definition: SVDBase.h:334
bool m_computeThinU
Definition: SVDBase.h:336
internal::variable_if_dynamic< Index, ColsAtCompileTime > m_cols
Definition: SVDBase.h:341
internal::variable_if_dynamic< Index, RowsAtCompileTime > m_rows
Definition: SVDBase.h:340
bool m_isInitialized
Definition: SVDBase.h:335
unsigned int m_computationOptions
Definition: SVDBase.h:338
MatrixVType m_matrixV
Definition: SVDBase.h:332
internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
Definition: SVDBase.h:158
internal::make_proper_matrix_type< Scalar, ColsAtCompileTime, MatrixVColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MatrixVMaxColsAtCompileTime >::type MatrixVType
Definition: SVDBase.h:156
Index m_nonzeroSingularValues
Definition: SVDBase.h:339
static constexpr bool ShouldComputeFullU
Definition: SVDBase.h:129
SingularValuesType m_singularValues
Definition: SVDBase.h:333
RealScalar m_prescribedThreshold
Definition: SVDBase.h:343
static constexpr bool ShouldComputeThinU
Definition: SVDBase.h:130
bool m_computeThinV
Definition: SVDBase.h:337
bool m_computeFullU
Definition: SVDBase.h:336
MatrixUType m_matrixU
Definition: SVDBase.h:331
bool m_isAllocated
Definition: SVDBase.h:335
@ Success
Definition: Constants.h:440

Member Function Documentation

◆ _check_compute_assertions()

template<typename Derived >
void Eigen::SVDBase< Derived >::_check_compute_assertions ( ) const
inlineprotected

◆ _check_solve_assertion()

template<typename Derived >
template<bool Transpose_, typename Rhs >
void Eigen::SVDBase< Derived >::_check_solve_assertion ( const Rhs &  b) const
inlineprotected
319  {
322  eigen_assert(computeU() && computeV() &&
323  "SVDBase::solve(): Both unitaries U and V are required to be computed (thin unitaries suffice).");
324  eigen_assert((Transpose_ ? cols() : rows()) == b.rows() &&
325  "SVDBase::solve(): invalid number of rows of the right hand side matrix b");
326  }
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
Definition: Macros.h:922
Scalar * b
Definition: benchVecAdd.cpp:17
Index cols() const
Definition: SVDBase.h:278
bool computeV() const
Definition: SVDBase.h:275
bool computeU() const
Definition: SVDBase.h:273
Index rows() const
Definition: SVDBase.h:277
void _check_compute_assertions() const
Definition: SVDBase.h:316

References Eigen::SVDBase< Derived >::_check_compute_assertions(), b, Eigen::SVDBase< Derived >::cols(), Eigen::SVDBase< Derived >::computeU(), Eigen::SVDBase< Derived >::computeV(), eigen_assert, EIGEN_ONLY_USED_FOR_DEBUG, and Eigen::SVDBase< Derived >::rows().

◆ _solve_impl()

template<typename Derived >
template<typename RhsType , typename DstType >
void Eigen::SVDBase< Derived >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const
372  {
373  // A = U S V^*
374  // So A^{-1} = V S^{-1} U^*
375 
376  Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime,
377  RhsType::MaxColsAtCompileTime>
378  tmp;
379  Index l_rank = rank();
380  tmp.noalias() = m_matrixU.leftCols(l_rank).adjoint() * rhs;
381  tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
382  dst = m_matrixV.leftCols(l_rank) * tmp;
383 }
SCALAR Scalar
Definition: bench_gemm.cpp:45
Index rank() const
Definition: SVDBase.h:217
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
const int Dynamic
Definition: Constants.h:25
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43

References Eigen::Dynamic, and tmp.

◆ _solve_impl_transposed()

template<typename Derived >
template<bool Conjugate, typename RhsType , typename DstType >
void Eigen::SVDBase< Derived >::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const
387  {
388  // A = U S V^*
389  // So A^{-*} = U S^{-1} V^*
390  // And A^{-T} = U_conj S^{-1} V^T
391  Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime,
392  RhsType::MaxColsAtCompileTime>
393  tmp;
394  Index l_rank = rank();
395 
396  tmp.noalias() = m_matrixV.leftCols(l_rank).transpose().template conjugateIf<Conjugate>() * rhs;
397  tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
398  dst = m_matrixU.template conjugateIf<!Conjugate>().leftCols(l_rank) * tmp;
399 }

References Eigen::Dynamic, and tmp.

◆ allocate()

template<typename Derived >
bool Eigen::SVDBase< Derived >::allocate ( Index  rows,
Index  cols,
unsigned int  computationOptions 
)
protected
403  {
404  eigen_assert(rows >= 0 && cols >= 0);
405 
406  if (m_isAllocated && rows == m_rows.value() && cols == m_cols.value() && computationOptions == m_computationOptions) {
407  return true;
408  }
409 
412  m_info = Success;
413  m_isInitialized = false;
414  m_isAllocated = true;
415  m_computationOptions = computationOptions;
420 
421  eigen_assert(!(m_computeFullU && m_computeThinU) && "SVDBase: you can't ask for both full and thin U");
422  eigen_assert(!(m_computeFullV && m_computeThinV) && "SVDBase: you can't ask for both full and thin V");
423 
426  if (RowsAtCompileTime == Dynamic)
428  if (ColsAtCompileTime == Dynamic)
430 
431  return false;
432 }
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void setValue(T v) const
Definition: XprHelper.h:163
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR T value()
Definition: XprHelper.h:161
constexpr bool should_svd_compute_full_v(int options)
Definition: SVDBase.h:39
constexpr bool should_svd_compute_thin_u(int options)
Definition: SVDBase.h:36
constexpr bool should_svd_compute_thin_v(int options)
Definition: SVDBase.h:38
constexpr bool should_svd_compute_full_u(int options)
Definition: SVDBase.h:37
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
Definition: MathFunctions.h:920

References cols, Eigen::Dynamic, eigen_assert, Eigen::numext::mini(), rows, Eigen::internal::should_svd_compute_full_u(), Eigen::internal::should_svd_compute_full_v(), Eigen::internal::should_svd_compute_thin_u(), Eigen::internal::should_svd_compute_thin_v(), and Eigen::Success.

Referenced by Eigen::JacobiSVD< MatrixType_, Options_ >::allocate().

◆ cols()

◆ computeU()

template<typename Derived >
bool Eigen::SVDBase< Derived >::computeU ( ) const
inline
Returns
true if U (full or thin) is asked for in this SVD decomposition
273 { return m_computeFullU || m_computeThinU; }

References Eigen::SVDBase< Derived >::m_computeFullU, and Eigen::SVDBase< Derived >::m_computeThinU.

Referenced by Eigen::SVDBase< Derived >::_check_solve_assertion(), and Eigen::SVDBase< Derived >::matrixU().

◆ computeV()

template<typename Derived >
bool Eigen::SVDBase< Derived >::computeV ( ) const
inline
Returns
true if V (full or thin) is asked for in this SVD decomposition
275 { return m_computeFullV || m_computeThinV; }

References Eigen::SVDBase< Derived >::m_computeFullV, and Eigen::SVDBase< Derived >::m_computeThinV.

Referenced by Eigen::SVDBase< Derived >::_check_solve_assertion(), and Eigen::SVDBase< Derived >::matrixV().

◆ derived() [1/2]

template<typename Derived >
Derived& Eigen::SVDBase< Derived >::derived ( )
inline
160 { return *static_cast<Derived*>(this); }

Referenced by Eigen::SVDBase< Derived >::setThreshold().

◆ derived() [2/2]

template<typename Derived >
const Derived& Eigen::SVDBase< Derived >::derived ( ) const
inline
161 { return *static_cast<const Derived*>(this); }

◆ diagSize()

template<typename Derived >
Index Eigen::SVDBase< Derived >::diagSize ( ) const
inline

◆ info()

template<typename Derived >
EIGEN_DEVICE_FUNC ComputationInfo Eigen::SVDBase< Derived >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful.
300  {
301  eigen_assert(m_isInitialized && "SVD is not initialized.");
302  return m_info;
303  }

References eigen_assert, Eigen::SVDBase< Derived >::m_info, and Eigen::SVDBase< Derived >::m_isInitialized.

◆ matrixU()

template<typename Derived >
const MatrixUType& Eigen::SVDBase< Derived >::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.

173  {
175  eigen_assert(computeU() && "This SVD decomposition didn't compute U. Did you ask for it?");
176  return m_matrixU;
177  }

References Eigen::SVDBase< Derived >::_check_compute_assertions(), Eigen::SVDBase< Derived >::computeU(), eigen_assert, and Eigen::SVDBase< Derived >::m_matrixU.

◆ matrixV()

template<typename Derived >
const MatrixVType& Eigen::SVDBase< Derived >::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.

189  {
191  eigen_assert(computeV() && "This SVD decomposition didn't compute V. Did you ask for it?");
192  return m_matrixV;
193  }

References Eigen::SVDBase< Derived >::_check_compute_assertions(), Eigen::SVDBase< Derived >::computeV(), eigen_assert, and Eigen::SVDBase< Derived >::m_matrixV.

◆ nonzeroSingularValues()

template<typename Derived >
Index Eigen::SVDBase< Derived >::nonzeroSingularValues ( ) const
inline
Returns
the number of singular values that are not exactly 0
206  {
209  }

References Eigen::SVDBase< Derived >::_check_compute_assertions(), and Eigen::SVDBase< Derived >::m_nonzeroSingularValues.

◆ rank()

template<typename Derived >
Index Eigen::SVDBase< Derived >::rank ( ) const
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&).
217  {
218  using std::abs;
220  if (m_singularValues.size() == 0) return 0;
221  RealScalar premultiplied_threshold =
222  numext::maxi<RealScalar>(m_singularValues.coeff(0) * threshold(), (std::numeric_limits<RealScalar>::min)());
224  while (i >= 0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
225  return i + 1;
226  }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
int i
Definition: BiCGSTAB_step_by_step.cpp:9
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
RealScalar threshold() const
Definition: SVDBase.h:265
#define min(a, b)
Definition: datatypes.h:22

References Eigen::SVDBase< Derived >::_check_compute_assertions(), abs(), i, Eigen::SVDBase< Derived >::m_nonzeroSingularValues, Eigen::SVDBase< Derived >::m_singularValues, min, and Eigen::SVDBase< Derived >::threshold().

◆ rows()

◆ setThreshold() [1/2]

template<typename Derived >
Derived& Eigen::SVDBase< Derived >::setThreshold ( const RealScalar threshold)
inline

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. This is not used for the SVD decomposition itself.

When it needs to get the threshold value, Eigen calls threshold(). The default is NumTraits<Scalar>::epsilon()

Parameters
thresholdThe new value to use as the threshold.

A singular value will be considered nonzero if its value is strictly greater than \( \vert singular value \vert \leqslant threshold \times \vert max singular value \vert \).

If you want to come back to the default behavior, call setThreshold(Default_t)

242  {
245  return derived();
246  }
Derived & derived()
Definition: SVDBase.h:160

References Eigen::SVDBase< Derived >::derived(), Eigen::SVDBase< Derived >::m_prescribedThreshold, Eigen::SVDBase< Derived >::m_usePrescribedThreshold, and Eigen::SVDBase< Derived >::threshold().

◆ setThreshold() [2/2]

template<typename Derived >
Derived& Eigen::SVDBase< Derived >::setThreshold ( Default_t  )
inline

Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold.

You should pass the special object Eigen::Default as parameter here.

svd.setThreshold(Eigen::Default);
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
@ Default
Definition: Constants.h:361

See the documentation of setThreshold(const RealScalar&).

256  {
257  m_usePrescribedThreshold = false;
258  return derived();
259  }

References Eigen::SVDBase< Derived >::derived(), and Eigen::SVDBase< Derived >::m_usePrescribedThreshold.

◆ singularValues()

template<typename Derived >
const SingularValuesType& Eigen::SVDBase< Derived >::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.

200  {
202  return m_singularValues;
203  }

References Eigen::SVDBase< Derived >::_check_compute_assertions(), and Eigen::SVDBase< Derived >::m_singularValues.

Referenced by check_generateRandomMatrixSvs(), compare_bdc_jacobi(), and Eigen::BDCSVD< MatrixType_, Options_ >::computeSVDofM().

◆ threshold()

template<typename Derived >
RealScalar Eigen::SVDBase< Derived >::threshold ( ) const
inline

Returns the threshold that will be used by certain methods such as rank().

See the documentation of setThreshold(const RealScalar&).

265  {
267  // this temporary is needed to workaround a MSVC issue
268  Index diagSize = (std::max<Index>)(1, m_diagSize);
270  }
NumTraits< typename MatrixType::Scalar >::Real RealScalar
Definition: SVDBase.h:126
Index diagSize() const
Definition: SVDBase.h:279
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43

References Eigen::SVDBase< Derived >::diagSize(), eigen_assert, oomph::SarahBL::epsilon, Eigen::SVDBase< Derived >::m_diagSize, Eigen::SVDBase< Derived >::m_isInitialized, Eigen::SVDBase< Derived >::m_prescribedThreshold, and Eigen::SVDBase< Derived >::m_usePrescribedThreshold.

Referenced by Eigen::SVDBase< Derived >::rank(), and Eigen::SVDBase< Derived >::setThreshold().

Friends And Related Function Documentation

◆ internal::solve_assertion

template<typename Derived >
template<typename Derived_ >
friend struct internal::solve_assertion
friend

Member Data Documentation

◆ m_cols

template<typename Derived >
internal::variable_if_dynamic<Index, ColsAtCompileTime> Eigen::SVDBase< Derived >::m_cols
protected

◆ m_computationOptions

template<typename Derived >
unsigned int Eigen::SVDBase< Derived >::m_computationOptions
protected

◆ m_computeFullU

template<typename Derived >
bool Eigen::SVDBase< Derived >::m_computeFullU
protected

◆ m_computeFullV

template<typename Derived >
bool Eigen::SVDBase< Derived >::m_computeFullV
protected

◆ m_computeThinU

template<typename Derived >
bool Eigen::SVDBase< Derived >::m_computeThinU
protected

◆ m_computeThinV

template<typename Derived >
bool Eigen::SVDBase< Derived >::m_computeThinV
protected

◆ m_diagSize

◆ m_info

template<typename Derived >
ComputationInfo Eigen::SVDBase< Derived >::m_info
protected

◆ m_isAllocated

template<typename Derived >
bool Eigen::SVDBase< Derived >::m_isAllocated
protected

◆ m_isInitialized

◆ m_matrixU

template<typename Derived >
MatrixUType Eigen::SVDBase< Derived >::m_matrixU
protected

◆ m_matrixV

template<typename Derived >
MatrixVType Eigen::SVDBase< Derived >::m_matrixV
protected

◆ m_nonzeroSingularValues

template<typename Derived >
Index Eigen::SVDBase< Derived >::m_nonzeroSingularValues
protected

◆ m_prescribedThreshold

template<typename Derived >
RealScalar Eigen::SVDBase< Derived >::m_prescribedThreshold
protected

◆ m_rows

template<typename Derived >
internal::variable_if_dynamic<Index, RowsAtCompileTime> Eigen::SVDBase< Derived >::m_rows
protected

◆ m_singularValues

template<typename Derived >
SingularValuesType Eigen::SVDBase< Derived >::m_singularValues
protected

◆ m_usePrescribedThreshold

template<typename Derived >
bool Eigen::SVDBase< Derived >::m_usePrescribedThreshold
protected

◆ ShouldComputeFullU

template<typename Derived >
constexpr bool Eigen::SVDBase< Derived >::ShouldComputeFullU = internal::traits<Derived>::ShouldComputeFullU
staticconstexpr

◆ ShouldComputeFullV

template<typename Derived >
constexpr bool Eigen::SVDBase< Derived >::ShouldComputeFullV = internal::traits<Derived>::ShouldComputeFullV
staticconstexpr

◆ ShouldComputeThinU

template<typename Derived >
constexpr bool Eigen::SVDBase< Derived >::ShouldComputeThinU = internal::traits<Derived>::ShouldComputeThinU
staticconstexpr

◆ ShouldComputeThinV

template<typename Derived >
constexpr bool Eigen::SVDBase< Derived >::ShouldComputeThinV = internal::traits<Derived>::ShouldComputeThinV
staticconstexpr

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