16 #ifndef EIGEN_SVDBASE_H
17 #define EIGEN_SVDBASE_H
40 template<
typename MatrixType,
int Options>
43 "SVDBase: Cannot request U or V using both static and runtime options, even if they match. "
44 "Requesting unitaries at runtime is DEPRECATED: "
45 "Prefer requesting unitaries statically, using the Options template parameter.");
48 "SVDBase: If thin U is requested at runtime, your matrix must have more rows than columns or a dynamic number of rows."
49 "Similarly, if thin V is requested at runtime, you matrix must have more columns than rows or a dynamic number of columns.");
50 (void)computationOptions;
55 template<
typename Derived>
struct traits<
SVDBase<Derived> >
60 typedef int StorageIndex;
64 template <
typename MatrixType,
int Options_>
65 struct svd_traits : traits<MatrixType> {
66 static constexpr
int Options = Options_;
72 DiagSizeAtCompileTime =
74 MaxDiagSizeAtCompileTime =
76 MatrixUColsAtCompileTime = ShouldComputeThinU ? DiagSizeAtCompileTime
77 : MatrixType::RowsAtCompileTime,
78 MatrixVColsAtCompileTime = ShouldComputeThinV ? DiagSizeAtCompileTime
79 : MatrixType::ColsAtCompileTime,
80 MatrixUMaxColsAtCompileTime = ShouldComputeThinU ? MaxDiagSizeAtCompileTime
81 : MatrixType::MaxRowsAtCompileTime,
82 MatrixVMaxColsAtCompileTime = ShouldComputeThinV ? MaxDiagSizeAtCompileTime
83 : MatrixType::MaxColsAtCompileTime
123 template<
typename Derived_>
124 friend struct internal::solve_assertion;
127 typedef typename MatrixType::Scalar
Scalar;
129 typedef typename Eigen::internal::traits<SVDBase>::StorageIndex
StorageIndex;
163 Derived&
derived() {
return *
static_cast<Derived*
>(
this); }
164 const Derived&
derived()
const {
return *
static_cast<const Derived*
>(
this); }
289 #ifdef EIGEN_PARSED_BY_DOXYGEN
299 template<
typename Rhs>
316 #ifndef EIGEN_PARSED_BY_DOXYGEN
317 template<
typename RhsType,
typename DstType>
318 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
320 template<
bool Conjugate,
typename RhsType,
typename DstType>
321 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const;
332 template<
bool Transpose_,
typename Rhs>
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");
378 #ifndef EIGEN_PARSED_BY_DOXYGEN
379 template<
typename Derived>
380 template<
typename RhsType,
typename DstType>
381 void SVDBase<Derived>::_solve_impl(
const RhsType &rhs, DstType &dst)
const
386 Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp;
387 Index l_rank = rank();
388 tmp.noalias() = m_matrixU.leftCols(l_rank).adjoint() * rhs;
389 tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
390 dst = m_matrixV.leftCols(l_rank) * tmp;
393 template<
typename Derived>
394 template<
bool Conjugate,
typename RhsType,
typename DstType>
395 void SVDBase<Derived>::_solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const
400 Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp;
401 Index l_rank = rank();
403 tmp.noalias() = m_matrixV.leftCols(l_rank).transpose().template conjugateIf<Conjugate>() * rhs;
404 tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
405 dst = m_matrixU.template conjugateIf<!Conjugate>().leftCols(l_rank) * tmp;
409 template <
typename Derived>
416 computationOptions == m_computationOptions)
424 m_isInitialized =
false;
425 m_isAllocated =
true;
426 m_computationOptions = computationOptions;
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");
435 m_diagSize = (
std::min)(m_rows, m_cols);
436 m_singularValues.resize(m_diagSize);
438 m_matrixU.resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0);
440 m_matrixV.resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0);
const AbsReturnType abs() const
#define EIGEN_DEVICE_FUNC
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
#define EIGEN_STATIC_ASSERT(X, MSG)
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Matrix< float, 1, Dynamic > MatrixType
Base class for all dense matrices, vectors, and expressions.
The matrix class, also used for vectors and row-vectors.
Base class of SVD algorithms.
void _check_solve_assertion(const Rhs &b) const
ComputationInfo info() const
Reports whether previous computation was successful.
static constexpr bool ShouldComputeFullV
Derived & setThreshold(const RealScalar &threshold)
bool m_usePrescribedThreshold
const MatrixVType & matrixV() const
Derived & setThreshold(Default_t)
static constexpr bool ShouldComputeThinV
const SingularValuesType & singularValues() const
NumTraits< typename MatrixType::Scalar >::Real RealScalar
unsigned int m_computationOptions
internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
RealScalar threshold() const
MatrixType::Scalar Scalar
Index m_nonzeroSingularValues
static constexpr bool ShouldComputeFullU
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
SingularValuesType m_singularValues
RealScalar m_prescribedThreshold
SVDBase()
Default Constructor.
static constexpr bool ShouldComputeThinU
bool allocate(Index rows, Index cols, unsigned int computationOptions)
@ MatrixVColsAtCompileTime
@ MatrixVMaxColsAtCompileTime
@ MaxDiagSizeAtCompileTime
@ MatrixUColsAtCompileTime
@ MatrixUMaxColsAtCompileTime
const Derived & derived() const
Eigen::internal::traits< SVDBase >::StorageIndex StorageIndex
internal::traits< Derived >::MatrixType MatrixType
void _check_compute_assertions() const
const MatrixUType & matrixU() const
Index nonzeroSingularValues() const
Pseudo expression representing a solving operation.
A base class for matrix decomposition and solvers.
@ HouseholderQRPreconditioner
@ ColPivHouseholderQRPreconditioner
@ FullPivHouseholderQRPreconditioner
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
constexpr bool should_svd_compute_full_v(int options)
constexpr bool should_svd_compute_thin_u(int options)
constexpr int get_qr_preconditioner(int options)
constexpr bool should_svd_compute_thin_v(int options)
constexpr int min_size_prefer_fixed(A a, B b)
constexpr int get_computation_options(int options)
void check_svd_options_assertions(unsigned int computationOptions, Index rows, Index cols)
constexpr bool should_svd_compute_full_u(int options)
constexpr int min_size_prefer_dynamic(A a, B b)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Eigen::Index Index
The interface type of indices.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.