![]() |
|
Base class of SVD algorithms. More...
#include <SVDBase.h>
Public Member Functions | |
Derived & | derived () |
const Derived & | derived () const |
const MatrixUType & | matrixU () const |
const MatrixVType & | matrixV () const |
const SingularValuesType & | singularValues () 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 |
![]() | |
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 |
![]() | |
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... | |
![]() | |
void | _check_solve_assertion (const Rhs &b) const |
Friends | |
template<typename Derived_ > | |
struct | internal::solve_assertion |
Base class of SVD algorithms.
Derived | the 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.
typedef internal::traits<Derived>::MatrixType Eigen::SVDBase< Derived >::MatrixType |
typedef NumTraits<typename MatrixType::Scalar>::Real Eigen::SVDBase< Derived >::RealScalar |
typedef MatrixType::Scalar Eigen::SVDBase< Derived >::Scalar |
typedef internal::plain_diag_type<MatrixType, RealScalar>::type Eigen::SVDBase< Derived >::SingularValuesType |
typedef Eigen::internal::traits<SVDBase>::StorageIndex Eigen::SVDBase< Derived >::StorageIndex |
anonymous enum |
|
inlineprotected |
Default Constructor.
Default constructor of SVDBase
|
inlineprotected |
References eigen_assert, and Eigen::SVDBase< Derived >::m_isInitialized.
Referenced by Eigen::SVDBase< Derived >::_check_solve_assertion(), Eigen::SVDBase< Derived >::matrixU(), Eigen::SVDBase< Derived >::matrixV(), Eigen::SVDBase< Derived >::nonzeroSingularValues(), Eigen::SVDBase< Derived >::rank(), and Eigen::SVDBase< Derived >::singularValues().
|
inlineprotected |
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().
void Eigen::SVDBase< Derived >::_solve_impl | ( | const RhsType & | rhs, |
DstType & | dst | ||
) | const |
References Eigen::Dynamic, and tmp.
void Eigen::SVDBase< Derived >::_solve_impl_transposed | ( | const RhsType & | rhs, |
DstType & | dst | ||
) | const |
References Eigen::Dynamic, and tmp.
|
protected |
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().
|
inline |
References Eigen::SVDBase< Derived >::m_cols, and Eigen::internal::variable_if_dynamic< T, Value >::value().
Referenced by gdb.printers._MatrixEntryIterator::__next__(), Eigen::SVDBase< Derived >::_check_solve_assertion(), gdb.printers.EigenMatrixPrinter::children(), gdb.printers.EigenSparseMatrixPrinter::children(), gdb.printers.EigenMatrixPrinter::to_string(), and gdb.printers.EigenSparseMatrixPrinter::to_string().
|
inline |
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().
|
inline |
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().
|
inline |
Referenced by Eigen::SVDBase< Derived >::setThreshold().
|
inline |
|
inline |
References Eigen::SVDBase< Derived >::m_diagSize, and Eigen::internal::variable_if_dynamic< T, Value >::value().
Referenced by Eigen::SVDBase< Derived >::threshold().
|
inline |
Reports whether previous computation was successful.
Success
if computation was successful. References eigen_assert, Eigen::SVDBase< Derived >::m_info, and Eigen::SVDBase< Derived >::m_isInitialized.
|
inline |
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.
References Eigen::SVDBase< Derived >::_check_compute_assertions(), Eigen::SVDBase< Derived >::computeU(), eigen_assert, and Eigen::SVDBase< Derived >::m_matrixU.
|
inline |
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.
References Eigen::SVDBase< Derived >::_check_compute_assertions(), Eigen::SVDBase< Derived >::computeV(), eigen_assert, and Eigen::SVDBase< Derived >::m_matrixV.
|
inline |
References Eigen::SVDBase< Derived >::_check_compute_assertions(), and Eigen::SVDBase< Derived >::m_nonzeroSingularValues.
|
inline |
*this
is the SVD.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().
|
inline |
References Eigen::SVDBase< Derived >::m_rows, and Eigen::internal::variable_if_dynamic< T, Value >::value().
Referenced by gdb.printers._MatrixEntryIterator::__next__(), Eigen::SVDBase< Derived >::_check_solve_assertion(), gdb.printers.EigenMatrixPrinter::children(), gdb.printers.EigenSparseMatrixPrinter::children(), gdb.printers.EigenMatrixPrinter::to_string(), and gdb.printers.EigenSparseMatrixPrinter::to_string().
|
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()
threshold | The 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)
References Eigen::SVDBase< Derived >::derived(), Eigen::SVDBase< Derived >::m_prescribedThreshold, Eigen::SVDBase< Derived >::m_usePrescribedThreshold, and Eigen::SVDBase< Derived >::threshold().
|
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.
See the documentation of setThreshold(const RealScalar&).
References Eigen::SVDBase< Derived >::derived(), and Eigen::SVDBase< Derived >::m_usePrescribedThreshold.
|
inline |
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.
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().
|
inline |
Returns the threshold that will be used by certain methods such as rank().
See the documentation of setThreshold(const RealScalar&).
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().
|
friend |
|
protected |
Referenced by Eigen::SVDBase< Derived >::cols().
|
protected |
|
protected |
Referenced by Eigen::SVDBase< Derived >::computeU().
|
protected |
Referenced by Eigen::SVDBase< Derived >::computeV().
|
protected |
Referenced by Eigen::SVDBase< Derived >::computeU().
|
protected |
Referenced by Eigen::SVDBase< Derived >::computeV().
|
protected |
Referenced by Eigen::SVDBase< Derived >::diagSize(), and Eigen::SVDBase< Derived >::threshold().
|
protected |
Referenced by Eigen::SVDBase< Derived >::info().
|
protected |
|
protected |
|
protected |
Referenced by Eigen::SVDBase< Derived >::matrixU().
|
protected |
Referenced by Eigen::SVDBase< Derived >::matrixV().
|
protected |
Referenced by Eigen::SVDBase< Derived >::nonzeroSingularValues(), and Eigen::SVDBase< Derived >::rank().
|
protected |
Referenced by Eigen::SVDBase< Derived >::setThreshold(), and Eigen::SVDBase< Derived >::threshold().
|
protected |
Referenced by Eigen::SVDBase< Derived >::rows().
|
protected |
Referenced by Eigen::SVDBase< Derived >::rank(), and Eigen::SVDBase< Derived >::singularValues().
|
protected |
Referenced by Eigen::SVDBase< Derived >::setThreshold(), and Eigen::SVDBase< Derived >::threshold().
|
staticconstexpr |
|
staticconstexpr |
|
staticconstexpr |
|
staticconstexpr |