Eigen::SVDBase< Derived > Class Template Reference

Base class of SVD algorithms. More...

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

Public Types

enum  {
  RowsAtCompileTime ,
  ColsAtCompileTime ,
  DiagSizeAtCompileTime ,
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime ,
  MaxDiagSizeAtCompileTime ,
  MatrixOptions ,
  MatrixUColsAtCompileTime ,
  MatrixVColsAtCompileTime ,
  MatrixUMaxColsAtCompileTime ,
  MatrixVMaxColsAtCompileTime
}
 
typedef Eigen::Index Index
 
typedef internal::traits< Derived >::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< SVDBase< Derived > >
enum  
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
typedef EigenBase< SVDBase< Derived > > Base
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const SVDBase< Derived > > ConstTransposeReturnType
 
typedef internal::traits< SVDBase< 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

Index cols () const
 
bool computeU () const
 
bool computeV () const
 
Derived & derived ()
 
const Derived & 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
 
Derived & setThreshold (const RealScalar &threshold)
 
Derived & setThreshold (Default_t)
 
const SingularValuesTypesingularValues () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
- Public Member Functions inherited from Eigen::SolverBase< SVDBase< Derived > >
const AdjointReturnType adjoint () const
 
SVDBase< Derived > & derived ()
 
const SVDBase< Derived > & derived () const
 
const Solve< SVDBase< 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
 

Static Public Attributes

static constexpr bool ShouldComputeFullU
 
static constexpr bool ShouldComputeFullV
 
static constexpr bool ShouldComputeThinU
 
static constexpr bool 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

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
 

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

Definition at line 118 of file SVDBase.h.

Member Typedef Documentation

◆ Index

template<typename Derived >
typedef Eigen::Index Eigen::SVDBase< Derived >::Index
Deprecated:
since Eigen 3.3

Definition at line 130 of file SVDBase.h.

◆ MatrixType

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

Definition at line 126 of file SVDBase.h.

◆ MatrixUType

template<typename Derived >
typedef internal::make_proper_matrix_type<Scalar, RowsAtCompileTime, MatrixUColsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MatrixUMaxColsAtCompileTime>::type Eigen::SVDBase< Derived >::MatrixUType

Definition at line 156 of file SVDBase.h.

◆ MatrixVType

template<typename Derived >
typedef internal::make_proper_matrix_type<Scalar, ColsAtCompileTime, MatrixVColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MatrixVMaxColsAtCompileTime>::type Eigen::SVDBase< Derived >::MatrixVType

Definition at line 159 of file SVDBase.h.

◆ RealScalar

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

Definition at line 128 of file SVDBase.h.

◆ Scalar

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

Definition at line 127 of file SVDBase.h.

◆ SingularValuesType

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

Definition at line 161 of file SVDBase.h.

◆ StorageIndex

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

Definition at line 129 of file SVDBase.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Derived >
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
DiagSizeAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
MaxDiagSizeAtCompileTime 
MatrixOptions 
MatrixUColsAtCompileTime 
MatrixVColsAtCompileTime 
MatrixUMaxColsAtCompileTime 
MatrixVMaxColsAtCompileTime 

Definition at line 137 of file SVDBase.h.

137  {
138  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
139  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
141  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
142  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
144  MatrixOptions = MatrixType::Options,
145  MatrixUColsAtCompileTime = internal::traits<Derived>::MatrixUColsAtCompileTime,
146  MatrixVColsAtCompileTime = internal::traits<Derived>::MatrixVColsAtCompileTime,
147  MatrixUMaxColsAtCompileTime = internal::traits<Derived>::MatrixUMaxColsAtCompileTime,
148  MatrixVMaxColsAtCompileTime = internal::traits<Derived>::MatrixVMaxColsAtCompileTime
149  };
@ MatrixVColsAtCompileTime
Definition: SVDBase.h:146
@ MatrixVMaxColsAtCompileTime
Definition: SVDBase.h:148
@ ColsAtCompileTime
Definition: SVDBase.h:139
@ DiagSizeAtCompileTime
Definition: SVDBase.h:140
@ MaxRowsAtCompileTime
Definition: SVDBase.h:141
@ MaxDiagSizeAtCompileTime
Definition: SVDBase.h:143
@ MatrixUColsAtCompileTime
Definition: SVDBase.h:145
@ MaxColsAtCompileTime
Definition: SVDBase.h:142
@ RowsAtCompileTime
Definition: SVDBase.h:138
@ MatrixUMaxColsAtCompileTime
Definition: SVDBase.h:147
constexpr int min_size_prefer_fixed(A a, B b)
Definition: Meta.h:553
constexpr int min_size_prefer_dynamic(A a, B b)
Definition: Meta.h:537

Constructor & Destructor Documentation

◆ SVDBase()

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

Default Constructor.

Default constructor of SVDBase

Definition at line 358 of file SVDBase.h.

359  : m_matrixU(MatrixUType()),
362  m_info(Success),
363  m_isInitialized(false),
364  m_isAllocated(false),
366  m_computeFullU(false),
367  m_computeThinU(false),
368  m_computeFullV(false),
369  m_computeThinV(false),
372  m_rows(-1),
373  m_cols(-1),
374  m_diagSize(0),
internal::make_proper_matrix_type< Scalar, RowsAtCompileTime, MatrixUColsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MatrixUMaxColsAtCompileTime >::type MatrixUType
Definition: SVDBase.h:156
bool m_usePrescribedThreshold
Definition: SVDBase.h:347
bool m_computeFullV
Definition: SVDBase.h:349
ComputationInfo m_info
Definition: SVDBase.h:346
bool m_computeThinU
Definition: SVDBase.h:348
bool m_isInitialized
Definition: SVDBase.h:347
Index m_cols
Definition: SVDBase.h:351
unsigned int m_computationOptions
Definition: SVDBase.h:350
MatrixVType m_matrixV
Definition: SVDBase.h:344
Index m_diagSize
Definition: SVDBase.h:351
internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
Definition: SVDBase.h:161
internal::make_proper_matrix_type< Scalar, ColsAtCompileTime, MatrixVColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MatrixVMaxColsAtCompileTime >::type MatrixVType
Definition: SVDBase.h:159
Index m_rows
Definition: SVDBase.h:351
Index m_nonzeroSingularValues
Definition: SVDBase.h:351
SingularValuesType m_singularValues
Definition: SVDBase.h:345
RealScalar m_prescribedThreshold
Definition: SVDBase.h:352
bool m_computeThinV
Definition: SVDBase.h:349
bool m_computeFullU
Definition: SVDBase.h:348
MatrixUType m_matrixU
Definition: SVDBase.h:343
bool m_isAllocated
Definition: SVDBase.h:347
@ Success
Definition: Constants.h:446

Member Function Documentation

◆ _check_compute_assertions()

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

Definition at line 328 of file SVDBase.h.

328  {
329  eigen_assert(m_isInitialized && "SVD is not initialized.");
330  }
#define eigen_assert(x)
Definition: Macros.h:902

◆ _check_solve_assertion()

template<typename Derived >
template<bool Transpose_, typename Rhs >
void Eigen::SVDBase< Derived >::_check_solve_assertion ( const Rhs &  b) const
inlineprotected

Definition at line 333 of file SVDBase.h.

333  {
336  eigen_assert(computeU() && computeV() && "SVDBase::solve(): Both unitaries U and V are required to be computed (thin unitaries suffice).");
337  eigen_assert((Transpose_?cols():rows())==b.rows() && "SVDBase::solve(): invalid number of rows of the right hand side matrix b");
338  }
Array< int, 3, 1 > b
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
Definition: Macros.h:914
Index cols() const
Definition: SVDBase.h:287
bool computeV() const
Definition: SVDBase.h:284
bool computeU() const
Definition: SVDBase.h:282
Index rows() const
Definition: SVDBase.h:286
void _check_compute_assertions() const
Definition: SVDBase.h:328

◆ allocate()

template<typename Derived >
bool Eigen::SVDBase< Derived >::allocate ( Index  rows,
Index  cols,
unsigned int  computationOptions 
)
protected

Definition at line 410 of file SVDBase.h.

410  {
411  eigen_assert(rows >= 0 && cols >= 0);
412 
413  if (m_isAllocated &&
414  rows == m_rows &&
415  cols == m_cols &&
416  computationOptions == m_computationOptions)
417  {
418  return true;
419  }
420 
421  m_rows = rows;
422  m_cols = cols;
423  m_info = Success;
424  m_isInitialized = false;
425  m_isAllocated = true;
426  m_computationOptions = computationOptions;
431 
432  eigen_assert(!(m_computeFullU && m_computeThinU) && "SVDBase: you can't ask for both full and thin U");
433  eigen_assert(!(m_computeFullV && m_computeThinV) && "SVDBase: you can't ask for both full and thin V");
434 
441 
442  return false;
443 }
constexpr void resize(Index rows, Index cols)
static constexpr bool ShouldComputeFullV
Definition: SVDBase.h:134
static constexpr bool ShouldComputeThinV
Definition: SVDBase.h:135
static constexpr bool ShouldComputeFullU
Definition: SVDBase.h:132
static constexpr bool ShouldComputeThinU
Definition: SVDBase.h:133
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
constexpr bool should_svd_compute_full_v(int options)
Definition: SVDBase.h:38
constexpr bool should_svd_compute_thin_u(int options)
Definition: SVDBase.h:35
constexpr bool should_svd_compute_thin_v(int options)
Definition: SVDBase.h:37
constexpr bool should_svd_compute_full_u(int options)
Definition: SVDBase.h:36
const int Dynamic
Definition: Constants.h:24

◆ cols()

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

Definition at line 287 of file SVDBase.h.

287 { return m_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

Definition at line 282 of file SVDBase.h.

282 { return m_computeFullU || m_computeThinU; }

◆ 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

Definition at line 284 of file SVDBase.h.

284 { return m_computeFullV || m_computeThinV; }

◆ derived() [1/2]

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

Definition at line 163 of file SVDBase.h.

163 { return *static_cast<Derived*>(this); }

◆ derived() [2/2]

template<typename Derived >
const Derived& Eigen::SVDBase< Derived >::derived ( ) const
inline

Definition at line 164 of file SVDBase.h.

164 { return *static_cast<const Derived*>(this); }

◆ info()

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

Reports whether previous computation was successful.

Returns
Success if computation was successful.

Definition at line 310 of file SVDBase.h.

311  {
312  eigen_assert(m_isInitialized && "SVD is not initialized.");
313  return m_info;
314  }

◆ 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.

Definition at line 175 of file SVDBase.h.

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

◆ 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.

Definition at line 191 of file SVDBase.h.

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

◆ nonzeroSingularValues()

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

Definition at line 210 of file SVDBase.h.

211  {
214  }

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

Definition at line 222 of file SVDBase.h.

223  {
224  using std::abs;
226  if(m_singularValues.size()==0) return 0;
227  RealScalar premultiplied_threshold = numext::maxi<RealScalar>(m_singularValues.coeff(0) * threshold(), (std::numeric_limits<RealScalar>::min)());
229  while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
230  return i+1;
231  }
const AbsReturnType abs() const
NumTraits< typename MatrixType::Scalar >::Real RealScalar
Definition: SVDBase.h:128
Eigen::Index Index
Definition: SVDBase.h:130
RealScalar threshold() const
Definition: SVDBase.h:272

◆ rows()

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

Definition at line 286 of file SVDBase.h.

286 { return m_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)

Definition at line 247 of file SVDBase.h.

248  {
251  return derived();
252  }
Derived & derived()
Definition: SVDBase.h:163

◆ 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:364

See the documentation of setThreshold(const RealScalar&).

Definition at line 262 of file SVDBase.h.

263  {
264  m_usePrescribedThreshold = false;
265  return derived();
266  }

◆ 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.

Definition at line 203 of file SVDBase.h.

204  {
206  return m_singularValues;
207  }

◆ solve()

template<typename Derived >
template<typename Rhs >
const Solve<Derived, Rhs> Eigen::SVDBase< Derived >::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 \).

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

Definition at line 272 of file SVDBase.h.

273  {
275  // this temporary is needed to workaround a MSVC issue
276  Index diagSize = (std::max<Index>)(1,m_diagSize);
278  : RealScalar(diagSize)*NumTraits<Scalar>::epsilon();
279  }

Member Data Documentation

◆ m_cols

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

Definition at line 351 of file SVDBase.h.

◆ m_computationOptions

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

Definition at line 350 of file SVDBase.h.

◆ m_computeFullU

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

Definition at line 348 of file SVDBase.h.

◆ m_computeFullV

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

Definition at line 349 of file SVDBase.h.

◆ m_computeThinU

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

Definition at line 348 of file SVDBase.h.

◆ m_computeThinV

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

Definition at line 349 of file SVDBase.h.

◆ m_diagSize

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

Definition at line 351 of file SVDBase.h.

◆ m_info

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

Definition at line 346 of file SVDBase.h.

◆ m_isAllocated

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

Definition at line 347 of file SVDBase.h.

◆ m_isInitialized

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

Definition at line 347 of file SVDBase.h.

◆ m_matrixU

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

Definition at line 343 of file SVDBase.h.

◆ m_matrixV

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

Definition at line 344 of file SVDBase.h.

◆ m_nonzeroSingularValues

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

Definition at line 351 of file SVDBase.h.

◆ m_prescribedThreshold

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

Definition at line 352 of file SVDBase.h.

◆ m_rows

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

Definition at line 351 of file SVDBase.h.

◆ m_singularValues

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

Definition at line 345 of file SVDBase.h.

◆ m_usePrescribedThreshold

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

Definition at line 347 of file SVDBase.h.

◆ ShouldComputeFullU

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

Definition at line 132 of file SVDBase.h.

◆ ShouldComputeFullV

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

Definition at line 134 of file SVDBase.h.

◆ ShouldComputeThinU

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

Definition at line 133 of file SVDBase.h.

◆ ShouldComputeThinV

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

Definition at line 135 of file SVDBase.h.


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