SVDBase.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009-2010 Benoit Jacob <jacob.benoit.1@gmail.com>
5 // Copyright (C) 2014 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // Copyright (C) 2013 Gauthier Brun <brun.gauthier@gmail.com>
8 // Copyright (C) 2013 Nicolas Carre <nicolas.carre@ensimag.fr>
9 // Copyright (C) 2013 Jean Ceccato <jean.ceccato@ensimag.fr>
10 // Copyright (C) 2013 Pierre Zoppitelli <pierre.zoppitelli@ensimag.fr>
11 //
12 // This Source Code Form is subject to the terms of the Mozilla
13 // Public License v. 2.0. If a copy of the MPL was not distributed
14 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
15 
16 #ifndef EIGEN_SVDBASE_H
17 #define EIGEN_SVDBASE_H
18 
19 // IWYU pragma: private
20 #include "./InternalHeaderCheck.h"
21 
22 namespace Eigen {
23 
24 namespace internal {
25 
30 };
31 
32 constexpr int get_qr_preconditioner(int options) { return options & QRPreconditionerBits; }
33 
34 constexpr int get_computation_options(int options) { return options & ComputationOptionsBits; }
35 
36 constexpr bool should_svd_compute_thin_u(int options) { return (options & ComputeThinU) != 0; }
37 constexpr bool should_svd_compute_full_u(int options) { return (options & ComputeFullU) != 0; }
38 constexpr bool should_svd_compute_thin_v(int options) { return (options & ComputeThinV) != 0; }
39 constexpr bool should_svd_compute_full_v(int options) { return (options & ComputeFullV) != 0; }
40 
41 template <typename MatrixType, int Options>
42 void check_svd_options_assertions(unsigned int computationOptions, Index rows, Index cols) {
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 "
51  "rows."
52  "Similarly, if thin V is requested at runtime, you matrix must have more columns than rows or a dynamic number "
53  "of columns.");
54  (void)computationOptions;
55  (void)rows;
56  (void)cols;
57 }
58 
59 template <typename Derived>
60 struct traits<SVDBase<Derived> > : traits<Derived> {
61  typedef MatrixXpr XprKind;
63  typedef int StorageIndex;
64  enum { Flags = 0 };
65 };
66 
67 template <typename MatrixType, int Options_>
68 struct svd_traits : traits<MatrixType> {
69  static constexpr int Options = Options_;
74  enum {
76  internal::min_size_prefer_dynamic(MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime),
78  internal::min_size_prefer_dynamic(MatrixType::MaxRowsAtCompileTime, MatrixType::MaxColsAtCompileTime),
82  MatrixVMaxColsAtCompileTime = ShouldComputeThinV ? MaxDiagSizeAtCompileTime : MatrixType::MaxColsAtCompileTime
83  };
84 };
85 } // namespace internal
86 
118 template <typename Derived>
119 class SVDBase : public SolverBase<SVDBase<Derived> > {
120  public:
121  template <typename Derived_>
123 
125  typedef typename MatrixType::Scalar Scalar;
128 
133 
134  enum {
135  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
136  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
138  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
139  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
146  };
147 
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")
150 
151  typedef
152  typename internal::make_proper_matrix_type<Scalar, RowsAtCompileTime, MatrixUColsAtCompileTime, MatrixOptions,
154  typedef
155  typename internal::make_proper_matrix_type<Scalar, ColsAtCompileTime, MatrixVColsAtCompileTime, MatrixOptions,
157 
158  typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
159 
160  Derived& derived() { return *static_cast<Derived*>(this); }
161  const Derived& derived() const { return *static_cast<const Derived*>(this); }
162 
173  const MatrixUType& matrixU() const {
175  eigen_assert(computeU() && "This SVD decomposition didn't compute U. Did you ask for it?");
176  return m_matrixU;
177  }
178 
189  const MatrixVType& matrixV() const {
191  eigen_assert(computeV() && "This SVD decomposition didn't compute V. Did you ask for it?");
192  return m_matrixV;
193  }
194 
202  return m_singularValues;
203  }
204 
209  }
210 
217  inline Index rank() const {
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  }
227 
242  Derived& setThreshold(const RealScalar& threshold) {
245  return derived();
246  }
247 
257  m_usePrescribedThreshold = false;
258  return derived();
259  }
260 
267  // this temporary is needed to workaround a MSVC issue
268  Index diagSize = (std::max<Index>)(1, m_diagSize);
270  }
271 
273  inline bool computeU() const { return m_computeFullU || m_computeThinU; }
275  inline bool computeV() const { return m_computeFullV || m_computeThinV; }
276 
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(); }
280 
281 #ifdef EIGEN_PARSED_BY_DOXYGEN
292  template <typename Rhs>
293  inline const Solve<Derived, Rhs> solve(const MatrixBase<Rhs>& b) const;
294 #endif
295 
301  eigen_assert(m_isInitialized && "SVD is not initialized.");
302  return m_info;
303  }
304 
305 #ifndef EIGEN_PARSED_BY_DOXYGEN
306  template <typename RhsType, typename DstType>
307  void _solve_impl(const RhsType& rhs, DstType& dst) const;
308 
309  template <bool Conjugate, typename RhsType, typename DstType>
310  void _solve_impl_transposed(const RhsType& rhs, DstType& dst) const;
311 #endif
312 
313  protected:
315 
316  void _check_compute_assertions() const { eigen_assert(m_isInitialized && "SVD is not initialized."); }
317 
318  template <bool Transpose_, typename Rhs>
319  void _check_solve_assertion(const Rhs& b) const {
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  }
327 
328  // return true if already allocated
329  bool allocate(Index rows, Index cols, unsigned int computationOptions);
330 
338  unsigned int m_computationOptions;
344 
350  : m_matrixU(MatrixUType()),
353  m_info(Success),
354  m_isInitialized(false),
355  m_isAllocated(false),
361  m_computationOptions(internal::traits<Derived>::Options),
367 };
368 
369 #ifndef EIGEN_PARSED_BY_DOXYGEN
370 template <typename Derived>
371 template <typename RhsType, typename DstType>
372 void SVDBase<Derived>::_solve_impl(const RhsType& rhs, DstType& dst) const {
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 }
384 
385 template <typename Derived>
386 template <bool Conjugate, typename RhsType, typename DstType>
387 void SVDBase<Derived>::_solve_impl_transposed(const RhsType& rhs, DstType& dst) const {
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 }
400 #endif
401 
402 template <typename Derived>
403 bool SVDBase<Derived>::allocate(Index rows, Index cols, unsigned int computationOptions) {
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 
410  m_rows.setValue(rows);
411  m_cols.setValue(cols);
412  m_info = Success;
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);
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 
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);
430 
431  return false;
432 }
433 
434 } // namespace Eigen
435 
436 #endif // EIGEN_SVDBASE_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
int i
Definition: BiCGSTAB_step_by_step.cpp:9
MatrixXcd V
Definition: EigenSolver_EigenSolver_MatrixType.cpp:15
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:892
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
Definition: Macros.h:922
#define eigen_assert(x)
Definition: Macros.h:910
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:74
int rows
Definition: Tutorial_commainit_02.cpp:1
int cols
Definition: Tutorial_commainit_02.cpp:1
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
Base class of SVD algorithms.
Definition: SVDBase.h:119
void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const
Definition: SVDBase.h:387
void _check_solve_assertion(const Rhs &b) const
Definition: SVDBase.h:319
Index cols() const
Definition: SVDBase.h:278
static constexpr bool ShouldComputeFullV
Definition: SVDBase.h:131
Derived & setThreshold(const RealScalar &threshold)
Definition: SVDBase.h:242
bool m_usePrescribedThreshold
Definition: SVDBase.h:335
bool m_computeFullV
Definition: SVDBase.h:337
const MatrixVType & matrixV() const
Definition: SVDBase.h:189
Derived & setThreshold(Default_t)
Definition: SVDBase.h:256
internal::variable_if_dynamic< Index, DiagSizeAtCompileTime > m_diagSize
Definition: SVDBase.h:342
static constexpr bool ShouldComputeThinV
Definition: SVDBase.h:132
Index rank() const
Definition: SVDBase.h:217
ComputationInfo m_info
Definition: SVDBase.h:334
bool m_computeThinU
Definition: SVDBase.h:336
const SingularValuesType & singularValues() const
Definition: SVDBase.h:200
internal::variable_if_dynamic< Index, ColsAtCompileTime > m_cols
Definition: SVDBase.h:341
EIGEN_DEVICE_FUNC ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SVDBase.h:300
NumTraits< typename MatrixType::Scalar >::Real RealScalar
Definition: SVDBase.h:126
bool computeV() const
Definition: SVDBase.h:275
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
bool computeU() const
Definition: SVDBase.h:273
internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
Definition: SVDBase.h:158
@ 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
RealScalar threshold() const
Definition: SVDBase.h:265
MatrixType::Scalar Scalar
Definition: SVDBase.h:125
Index m_nonzeroSingularValues
Definition: SVDBase.h:339
static constexpr bool ShouldComputeFullU
Definition: SVDBase.h:129
Index rows() const
Definition: SVDBase.h:277
Derived & derived()
Definition: SVDBase.h:160
SingularValuesType m_singularValues
Definition: SVDBase.h:333
RealScalar m_prescribedThreshold
Definition: SVDBase.h:343
SVDBase()
Default Constructor.
Definition: SVDBase.h:349
void _solve_impl(const RhsType &rhs, DstType &dst) const
Definition: SVDBase.h:372
static constexpr bool ShouldComputeThinU
Definition: SVDBase.h:130
bool m_computeThinV
Definition: SVDBase.h:337
bool m_computeFullU
Definition: SVDBase.h:336
bool allocate(Index rows, Index cols, unsigned int computationOptions)
Definition: SVDBase.h:403
Index diagSize() const
Definition: SVDBase.h:279
const Derived & derived() const
Definition: SVDBase.h:161
Eigen::internal::traits< SVDBase >::StorageIndex StorageIndex
Definition: SVDBase.h:127
internal::traits< Derived >::MatrixType MatrixType
Definition: SVDBase.h:124
MatrixUType m_matrixU
Definition: SVDBase.h:331
void _check_compute_assertions() const
Definition: SVDBase.h:316
bool m_isAllocated
Definition: SVDBase.h:335
const MatrixUType & matrixU() const
Definition: SVDBase.h:173
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
internal::traits< Derived >::Scalar Scalar
Definition: SolverBase.h:75
const Solve< SVDBase< Derived >, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SolverBase.h:106
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR T value()
Definition: XprHelper.h:161
#define min(a, b)
Definition: datatypes.h:22
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
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
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
@ Rhs
Definition: TensorContractionMapper.h:20
constexpr int get_qr_preconditioner(int options)
Definition: SVDBase.h:32
constexpr bool should_svd_compute_thin_v(int options)
Definition: SVDBase.h:38
constexpr int min_size_prefer_fixed(A a, B b)
Definition: Meta.h:683
OptionsMasks
Definition: SVDBase.h:26
@ QRPreconditionerBits
Definition: SVDBase.h:27
@ ComputationOptionsBits
Definition: SVDBase.h:29
constexpr int get_computation_options(int options)
Definition: SVDBase.h:34
void check_svd_options_assertions(unsigned int computationOptions, Index rows, Index cols)
Definition: SVDBase.h:42
constexpr bool should_svd_compute_full_u(int options)
Definition: SVDBase.h:37
constexpr int min_size_prefer_dynamic(A a, B b)
Definition: Meta.h:668
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
Definition: MathFunctions.h:920
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
Default_t
Definition: Constants.h:361
const int Dynamic
Definition: Constants.h:25
Extend namespace for flags.
Definition: fsi_chan_precond_driver.cc:56
type
Definition: compute_granudrum_aor.py:141
Definition: Eigen_Colamd.h:49
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
Definition: Constants.h:534
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: Constants.h:525
Definition: SolverBase.h:21
Definition: SVDBase.h:68
@ MatrixUMaxColsAtCompileTime
Definition: SVDBase.h:81
@ MatrixVColsAtCompileTime
Definition: SVDBase.h:80
@ MaxDiagSizeAtCompileTime
Definition: SVDBase.h:77
@ MatrixVMaxColsAtCompileTime
Definition: SVDBase.h:82
@ MatrixUColsAtCompileTime
Definition: SVDBase.h:79
@ DiagSizeAtCompileTime
Definition: SVDBase.h:75
static constexpr int Options
Definition: SVDBase.h:69
static constexpr bool ShouldComputeThinV
Definition: SVDBase.h:73
static constexpr bool ShouldComputeFullU
Definition: SVDBase.h:70
static constexpr bool ShouldComputeThinU
Definition: SVDBase.h:71
static constexpr bool ShouldComputeFullV
Definition: SVDBase.h:72
SolverStorage StorageKind
Definition: SVDBase.h:62
MatrixXpr XprKind
Definition: SVDBase.h:61
int StorageIndex
Definition: SVDBase.h:63
Definition: ForwardDeclarations.h:21