16#ifndef EIGEN_SVDBASE_H
17#define EIGEN_SVDBASE_H
20#include "./InternalHeaderCheck.h"
32constexpr int get_qr_preconditioner(
int options) {
return options & QRPreconditionerBits; }
34constexpr int get_computation_options(
int options) {
return options & ComputationOptionsBits; }
36constexpr bool should_svd_compute_thin_u(
int options) {
return (options &
ComputeThinU) != 0; }
37constexpr bool should_svd_compute_full_u(
int options) {
return (options &
ComputeFullU) != 0; }
38constexpr bool should_svd_compute_thin_v(
int options) {
return (options &
ComputeThinV) != 0; }
39constexpr bool should_svd_compute_full_v(
int options) {
return (options &
ComputeFullV) != 0; }
41template <
typename MatrixType,
int Options>
42void check_svd_options_assertions(
unsigned int computationOptions, Index rows, Index cols) {
43 EIGEN_STATIC_ASSERT((Options & ComputationOptionsBits) == 0,
44 "SVDBase: Cannot request U or V using both static and runtime options, even if they match. "
45 "Requesting unitaries at runtime is DEPRECATED: "
46 "Prefer requesting unitaries statically, using the Options template parameter.");
48 !(should_svd_compute_thin_u(computationOptions) && cols < rows && MatrixType::RowsAtCompileTime !=
Dynamic) &&
49 !(should_svd_compute_thin_v(computationOptions) && rows < cols && MatrixType::ColsAtCompileTime !=
Dynamic) &&
50 "SVDBase: If thin U is requested at runtime, your matrix must have more rows than columns or a dynamic number of "
52 "Similarly, if thin V is requested at runtime, you matrix must have more columns than rows or a dynamic number "
54 (void)computationOptions;
59template <
typename Derived>
60struct traits<SVDBase<Derived> > : traits<Derived> {
61 typedef MatrixXpr XprKind;
62 typedef SolverStorage StorageKind;
63 typedef int StorageIndex;
67template <
typename MatrixType,
int Options_>
68struct svd_traits : traits<MatrixType> {
69 static constexpr int Options = Options_;
70 static constexpr bool ShouldComputeFullU = internal::should_svd_compute_full_u(Options);
71 static constexpr bool ShouldComputeThinU = internal::should_svd_compute_thin_u(Options);
72 static constexpr bool ShouldComputeFullV = internal::should_svd_compute_full_v(Options);
73 static constexpr bool ShouldComputeThinV = internal::should_svd_compute_thin_v(Options);
75 DiagSizeAtCompileTime =
76 internal::min_size_prefer_dynamic(MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime),
77 MaxDiagSizeAtCompileTime =
78 internal::min_size_prefer_dynamic(MatrixType::MaxRowsAtCompileTime, MatrixType::MaxColsAtCompileTime),
79 MatrixUColsAtCompileTime = ShouldComputeThinU ? DiagSizeAtCompileTime : MatrixType::RowsAtCompileTime,
80 MatrixVColsAtCompileTime = ShouldComputeThinV ? DiagSizeAtCompileTime : MatrixType::ColsAtCompileTime,
81 MatrixUMaxColsAtCompileTime = ShouldComputeThinU ? MaxDiagSizeAtCompileTime : MatrixType::MaxRowsAtCompileTime,
82 MatrixVMaxColsAtCompileTime = ShouldComputeThinV ? MaxDiagSizeAtCompileTime : MatrixType::MaxColsAtCompileTime
118template <
typename Derived>
121 template <
typename Derived_>
122 friend struct internal::solve_assertion;
124 typedef typename internal::traits<Derived>::MatrixType MatrixType;
125 typedef typename MatrixType::Scalar Scalar;
126 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
127 typedef typename Eigen::internal::traits<SVDBase>::StorageIndex StorageIndex;
129 static constexpr bool ShouldComputeFullU = internal::traits<Derived>::ShouldComputeFullU;
130 static constexpr bool ShouldComputeThinU = internal::traits<Derived>::ShouldComputeThinU;
131 static constexpr bool ShouldComputeFullV = internal::traits<Derived>::ShouldComputeFullV;
132 static constexpr bool ShouldComputeThinV = internal::traits<Derived>::ShouldComputeThinV;
135 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
136 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
137 DiagSizeAtCompileTime = internal::min_size_prefer_dynamic(RowsAtCompileTime, ColsAtCompileTime),
138 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
139 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
140 MaxDiagSizeAtCompileTime = internal::min_size_prefer_fixed(MaxRowsAtCompileTime, 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
148 EIGEN_STATIC_ASSERT(!(ShouldComputeFullU && ShouldComputeThinU),
"SVDBase: Cannot request both full and thin U")
149 EIGEN_STATIC_ASSERT(!(ShouldComputeFullV && ShouldComputeThinV),
"SVDBase: Cannot request both full and thin V")
152 typename internal::make_proper_matrix_type<Scalar, RowsAtCompileTime, MatrixUColsAtCompileTime, MatrixOptions,
153 MaxRowsAtCompileTime, MatrixUMaxColsAtCompileTime>::type
MatrixUType;
155 typename internal::make_proper_matrix_type<Scalar, ColsAtCompileTime, MatrixVColsAtCompileTime, MatrixOptions,
156 MaxColsAtCompileTime, MatrixVMaxColsAtCompileTime>::type
MatrixVType;
158 typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
160 Derived& derived() {
return *
static_cast<Derived*
>(
this); }
161 const Derived& derived()
const {
return *
static_cast<const Derived*
>(
this); }
174 _check_compute_assertions();
175 eigen_assert(
computeU() &&
"This SVD decomposition didn't compute U. Did you ask for it?");
190 _check_compute_assertions();
191 eigen_assert(
computeV() &&
"This SVD decomposition didn't compute V. Did you ask for it?");
201 _check_compute_assertions();
202 return m_singularValues;
207 _check_compute_assertions();
208 return m_nonzeroSingularValues;
219 _check_compute_assertions();
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)());
223 Index i = m_nonzeroSingularValues - 1;
224 while (i >= 0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
243 m_usePrescribedThreshold =
true;
257 m_usePrescribedThreshold =
false;
266 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
268 Index diagSize = (std::max<Index>)(1, m_diagSize);
273 inline bool computeU()
const {
return m_computeFullU || m_computeThinU; }
275 inline bool computeV()
const {
return m_computeFullV || m_computeThinV; }
277 inline Index rows()
const {
return m_rows.value(); }
278 inline Index cols()
const {
return m_cols.value(); }
279 inline Index diagSize()
const {
return m_diagSize.value(); }
281#ifdef EIGEN_PARSED_BY_DOXYGEN
292 template <
typename Rhs>
301 eigen_assert(m_isInitialized &&
"SVD is not initialized.");
305#ifndef EIGEN_PARSED_BY_DOXYGEN
306 template <
typename RhsType,
typename DstType>
307 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
309 template <
bool Conjugate,
typename RhsType,
typename DstType>
310 void _solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const;
314 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
316 void _check_compute_assertions()
const { eigen_assert(m_isInitialized &&
"SVD is not initialized."); }
318 template <
bool Transpose_,
typename Rhs>
319 void _check_solve_assertion(
const Rhs& b)
const {
320 EIGEN_ONLY_USED_FOR_DEBUG(b);
321 _check_compute_assertions();
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");
329 bool allocate(Index rows, Index cols,
unsigned int computationOptions);
331 MatrixUType m_matrixU;
332 MatrixVType m_matrixV;
333 SingularValuesType m_singularValues;
334 ComputationInfo m_info;
335 bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
336 bool m_computeFullU, m_computeThinU;
337 bool m_computeFullV, m_computeThinV;
338 unsigned int m_computationOptions;
339 Index m_nonzeroSingularValues;
340 internal::variable_if_dynamic<Index, RowsAtCompileTime> m_rows;
341 internal::variable_if_dynamic<Index, ColsAtCompileTime> m_cols;
342 internal::variable_if_dynamic<Index, DiagSizeAtCompileTime> m_diagSize;
343 RealScalar m_prescribedThreshold;
352 m_singularValues(SingularValuesType()),
354 m_isInitialized(false),
355 m_isAllocated(false),
356 m_usePrescribedThreshold(false),
357 m_computeFullU(ShouldComputeFullU),
358 m_computeThinU(ShouldComputeThinU),
359 m_computeFullV(ShouldComputeFullV),
360 m_computeThinV(ShouldComputeThinV),
361 m_computationOptions(internal::traits<Derived>::Options),
362 m_nonzeroSingularValues(0),
363 m_rows(RowsAtCompileTime),
364 m_cols(ColsAtCompileTime),
365 m_diagSize(DiagSizeAtCompileTime),
366 m_prescribedThreshold(0) {}
369#ifndef EIGEN_PARSED_BY_DOXYGEN
370template <
typename Derived>
371template <
typename RhsType,
typename DstType>
372void SVDBase<Derived>::_solve_impl(
const RhsType& rhs, DstType& dst)
const {
376 Matrix<
typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime,
377 RhsType::MaxColsAtCompileTime>
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;
385template <
typename Derived>
386template <
bool Conjugate,
typename RhsType,
typename DstType>
387void SVDBase<Derived>::_solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const {
391 Matrix<
typename RhsType::Scalar,
Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime,
392 RhsType::MaxColsAtCompileTime>
394 Index l_rank = rank();
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;
402template <
typename Derived>
403bool SVDBase<Derived>::allocate(Index rows, Index cols,
unsigned int computationOptions) {
404 eigen_assert(rows >= 0 && cols >= 0);
406 if (m_isAllocated && rows == m_rows.value() && cols == m_cols.value() && computationOptions == m_computationOptions) {
410 m_rows.setValue(rows);
411 m_cols.setValue(cols);
413 m_isInitialized =
false;
414 m_isAllocated =
true;
415 m_computationOptions = computationOptions;
416 m_computeFullU = ShouldComputeFullU || internal::should_svd_compute_full_u(computationOptions);
417 m_computeThinU = ShouldComputeThinU || internal::should_svd_compute_thin_u(computationOptions);
418 m_computeFullV = ShouldComputeFullV || internal::should_svd_compute_full_v(computationOptions);
419 m_computeThinV = ShouldComputeThinV || internal::should_svd_compute_thin_v(computationOptions);
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");
424 m_diagSize.setValue(numext::mini(m_rows.value(), m_cols.value()));
425 m_singularValues.resize(m_diagSize.value());
426 if (RowsAtCompileTime == Dynamic)
427 m_matrixU.resize(m_rows.value(), m_computeFullU ? m_rows.value() : m_computeThinU ? m_diagSize.value() : 0);
428 if (ColsAtCompileTime == Dynamic)
429 m_matrixV.resize(m_cols.value(), m_computeFullV ? m_cols.value() : m_computeThinV ? m_diagSize.value() : 0);
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
Base class of SVD algorithms.
Definition SVDBase.h:119
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SVDBase.h:300
Derived & setThreshold(const RealScalar &threshold)
Definition SVDBase.h:242
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Index rank() const
Definition SVDBase.h:217
bool computeV() const
Definition SVDBase.h:275
bool computeU() const
Definition SVDBase.h:273
Derived & setThreshold(Default_t)
Definition SVDBase.h:256
RealScalar threshold() const
Definition SVDBase.h:265
SVDBase()
Default Constructor.
Definition SVDBase.h:349
const SingularValuesType & singularValues() const
Definition SVDBase.h:200
const MatrixUType & matrixU() const
Definition SVDBase.h:173
const MatrixVType & matrixV() const
Definition SVDBase.h:189
Index nonzeroSingularValues() const
Definition SVDBase.h:206
Pseudo expression representing a solving operation.
Definition Solve.h:62
A base class for matrix decomposition and solvers.
Definition SolverBase.h:72
ComputationInfo
Definition Constants.h:438
@ NoQRPreconditioner
Definition Constants.h:423
@ HouseholderQRPreconditioner
Definition Constants.h:425
@ ColPivHouseholderQRPreconditioner
Definition Constants.h:421
@ FullPivHouseholderQRPreconditioner
Definition Constants.h:427
@ Success
Definition Constants.h:440
@ ComputeFullV
Definition Constants.h:393
@ ComputeThinV
Definition Constants.h:395
@ ComputeFullU
Definition Constants.h:389
@ ComputeThinU
Definition Constants.h:391
Namespace containing all symbols from the Eigen library.
Definition Core:137
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83
const int Dynamic
Definition Constants.h:25
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523