Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
SVDBase.h
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 #include "./InternalHeaderCheck.h"
20 
21 namespace Eigen {
22 
23 namespace internal {
24 
25 enum OptionsMasks {
28  ComputationOptionsBits = ComputeThinU | ComputeFullU | ComputeThinV | ComputeFullV
29 };
30 
31 constexpr int get_qr_preconditioner(int options) { return options & QRPreconditionerBits; }
32 
33 constexpr int get_computation_options(int options) { return options & ComputationOptionsBits; }
34 
35 constexpr bool should_svd_compute_thin_u(int options) { return (options & ComputeThinU) != 0; }
36 constexpr bool should_svd_compute_full_u(int options) { return (options & ComputeFullU) != 0; }
37 constexpr bool should_svd_compute_thin_v(int options) { return (options & ComputeThinV) != 0; }
38 constexpr bool should_svd_compute_full_v(int options) { return (options & ComputeFullV) != 0; }
39 
40 template<typename MatrixType, int Options>
41 void check_svd_options_assertions(unsigned int computationOptions, Index rows, Index cols) {
42  EIGEN_STATIC_ASSERT((Options & ComputationOptionsBits) == 0,
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.");
46  eigen_assert(!(should_svd_compute_thin_u(computationOptions) && cols < rows && MatrixType::RowsAtCompileTime != Dynamic) &&
47  !(should_svd_compute_thin_v(computationOptions) && rows < cols && MatrixType::ColsAtCompileTime != Dynamic) &&
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;
51  (void)rows;
52  (void)cols;
53 }
54 
55 template<typename Derived> struct traits<SVDBase<Derived> >
56  : traits<Derived>
57 {
58  typedef MatrixXpr XprKind;
59  typedef SolverStorage StorageKind;
60  typedef int StorageIndex;
61  enum { Flags = 0 };
62 };
63 
64 template <typename MatrixType, int Options_>
65 struct svd_traits : traits<MatrixType> {
66  static constexpr int Options = Options_;
67  static constexpr bool ShouldComputeFullU = internal::should_svd_compute_full_u(Options);
68  static constexpr bool ShouldComputeThinU = internal::should_svd_compute_thin_u(Options);
69  static constexpr bool ShouldComputeFullV = internal::should_svd_compute_full_v(Options);
70  static constexpr bool ShouldComputeThinV = internal::should_svd_compute_thin_v(Options);
71  enum {
72  DiagSizeAtCompileTime =
73  internal::min_size_prefer_dynamic(MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime),
74  MaxDiagSizeAtCompileTime =
75  internal::min_size_prefer_dynamic(MatrixType::MaxRowsAtCompileTime, MatrixType::MaxColsAtCompileTime),
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
84  };
85 };
86 }
87 
118 template<typename Derived> class SVDBase
119  : public SolverBase<SVDBase<Derived> >
120 {
121 public:
122 
123  template<typename Derived_>
124  friend struct internal::solve_assertion;
125 
126  typedef typename internal::traits<Derived>::MatrixType MatrixType;
127  typedef typename MatrixType::Scalar Scalar;
128  typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
129  typedef typename Eigen::internal::traits<SVDBase>::StorageIndex StorageIndex;
130  typedef Eigen::Index Index;
131 
132  static constexpr bool ShouldComputeFullU = internal::traits<Derived>::ShouldComputeFullU;
133  static constexpr bool ShouldComputeThinU = internal::traits<Derived>::ShouldComputeThinU;
134  static constexpr bool ShouldComputeFullV = internal::traits<Derived>::ShouldComputeFullV;
135  static constexpr bool ShouldComputeThinV = internal::traits<Derived>::ShouldComputeThinV;
136 
137  enum {
138  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
139  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
140  DiagSizeAtCompileTime = internal::min_size_prefer_dynamic(RowsAtCompileTime, ColsAtCompileTime),
141  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
142  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
143  MaxDiagSizeAtCompileTime = internal::min_size_prefer_fixed(MaxRowsAtCompileTime, 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  };
150 
151  EIGEN_STATIC_ASSERT(!(ShouldComputeFullU && ShouldComputeThinU), "SVDBase: Cannot request both full and thin U")
152  EIGEN_STATIC_ASSERT(!(ShouldComputeFullV && ShouldComputeThinV), "SVDBase: Cannot request both full and thin V")
153 
154  typedef
155  typename internal::make_proper_matrix_type<Scalar, RowsAtCompileTime, MatrixUColsAtCompileTime, MatrixOptions,
156  MaxRowsAtCompileTime, MatrixUMaxColsAtCompileTime>::type MatrixUType;
157  typedef
158  typename internal::make_proper_matrix_type<Scalar, ColsAtCompileTime, MatrixVColsAtCompileTime, MatrixOptions,
159  MaxColsAtCompileTime, MatrixVMaxColsAtCompileTime>::type MatrixVType;
160 
161  typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
162 
163  Derived& derived() { return *static_cast<Derived*>(this); }
164  const Derived& derived() const { return *static_cast<const Derived*>(this); }
165 
175  const MatrixUType& matrixU() const
176  {
177  _check_compute_assertions();
178  eigen_assert(computeU() && "This SVD decomposition didn't compute U. Did you ask for it?");
179  return m_matrixU;
180  }
181 
191  const MatrixVType& matrixV() const
192  {
193  _check_compute_assertions();
194  eigen_assert(computeV() && "This SVD decomposition didn't compute V. Did you ask for it?");
195  return m_matrixV;
196  }
197 
203  const SingularValuesType& singularValues() const
204  {
205  _check_compute_assertions();
206  return m_singularValues;
207  }
208 
211  {
212  _check_compute_assertions();
213  return m_nonzeroSingularValues;
214  }
215 
222  inline Index rank() const
223  {
224  using std::abs;
225  _check_compute_assertions();
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)());
228  Index i = m_nonzeroSingularValues-1;
229  while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
230  return i+1;
231  }
232 
247  Derived& setThreshold(const RealScalar& threshold)
248  {
249  m_usePrescribedThreshold = true;
250  m_prescribedThreshold = threshold;
251  return derived();
252  }
253 
262  Derived& setThreshold(Default_t)
263  {
264  m_usePrescribedThreshold = false;
265  return derived();
266  }
267 
272  RealScalar threshold() const
273  {
274  eigen_assert(m_isInitialized || m_usePrescribedThreshold);
275  // this temporary is needed to workaround a MSVC issue
276  Index diagSize = (std::max<Index>)(1,m_diagSize);
277  return m_usePrescribedThreshold ? m_prescribedThreshold
278  : RealScalar(diagSize)*NumTraits<Scalar>::epsilon();
279  }
280 
282  inline bool computeU() const { return m_computeFullU || m_computeThinU; }
284  inline bool computeV() const { return m_computeFullV || m_computeThinV; }
285 
286  inline Index rows() const { return m_rows; }
287  inline Index cols() const { return m_cols; }
288 
289  #ifdef EIGEN_PARSED_BY_DOXYGEN
299  template<typename Rhs>
300  inline const Solve<Derived, Rhs>
301  solve(const MatrixBase<Rhs>& b) const;
302  #endif
303 
304 
309  EIGEN_DEVICE_FUNC
311  {
312  eigen_assert(m_isInitialized && "SVD is not initialized.");
313  return m_info;
314  }
315 
316  #ifndef EIGEN_PARSED_BY_DOXYGEN
317  template<typename RhsType, typename DstType>
318  void _solve_impl(const RhsType &rhs, DstType &dst) const;
319 
320  template<bool Conjugate, typename RhsType, typename DstType>
321  void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
322  #endif
323 
324 protected:
325 
326  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
327 
328  void _check_compute_assertions() const {
329  eigen_assert(m_isInitialized && "SVD is not initialized.");
330  }
331 
332  template<bool Transpose_, typename Rhs>
333  void _check_solve_assertion(const Rhs& b) const {
334  EIGEN_ONLY_USED_FOR_DEBUG(b);
335  _check_compute_assertions();
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  }
339 
340  // return true if already allocated
341  bool allocate(Index rows, Index cols, unsigned int computationOptions);
342 
343  MatrixUType m_matrixU;
344  MatrixVType m_matrixV;
345  SingularValuesType m_singularValues;
346  ComputationInfo m_info;
347  bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
348  bool m_computeFullU, m_computeThinU;
349  bool m_computeFullV, m_computeThinV;
350  unsigned int m_computationOptions;
351  Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
352  RealScalar m_prescribedThreshold;
353 
359  : m_info(Success),
360  m_isInitialized(false),
361  m_isAllocated(false),
362  m_usePrescribedThreshold(false),
363  m_computeFullU(false),
364  m_computeThinU(false),
365  m_computeFullV(false),
366  m_computeThinV(false),
367  m_computationOptions(0),
368  m_rows(-1),
369  m_cols(-1),
370  m_diagSize(0) {}
371 };
372 
373 #ifndef EIGEN_PARSED_BY_DOXYGEN
374 template<typename Derived>
375 template<typename RhsType, typename DstType>
376 void SVDBase<Derived>::_solve_impl(const RhsType &rhs, DstType &dst) const
377 {
378  // A = U S V^*
379  // So A^{-1} = V S^{-1} U^*
380 
381  Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp;
382  Index l_rank = rank();
383  tmp.noalias() = m_matrixU.leftCols(l_rank).adjoint() * rhs;
384  tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
385  dst = m_matrixV.leftCols(l_rank) * tmp;
386 }
387 
388 template<typename Derived>
389 template<bool Conjugate, typename RhsType, typename DstType>
390 void SVDBase<Derived>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
391 {
392  // A = U S V^*
393  // So A^{-*} = U S^{-1} V^*
394  // And A^{-T} = U_conj S^{-1} V^T
395  Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp;
396  Index l_rank = rank();
397 
398  tmp.noalias() = m_matrixV.leftCols(l_rank).transpose().template conjugateIf<Conjugate>() * rhs;
399  tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
400  dst = m_matrixU.template conjugateIf<!Conjugate>().leftCols(l_rank) * tmp;
401 }
402 #endif
403 
404 template <typename Derived>
405 bool SVDBase<Derived>::allocate(Index rows, Index cols, unsigned int computationOptions) {
406  eigen_assert(rows >= 0 && cols >= 0);
407 
408  if (m_isAllocated &&
409  rows == m_rows &&
410  cols == m_cols &&
411  computationOptions == m_computationOptions)
412  {
413  return true;
414  }
415 
416  m_rows = rows;
417  m_cols = cols;
418  m_info = Success;
419  m_isInitialized = false;
420  m_isAllocated = true;
421  m_computationOptions = computationOptions;
422  m_computeFullU = ShouldComputeFullU || internal::should_svd_compute_full_u(computationOptions);
423  m_computeThinU = ShouldComputeThinU || internal::should_svd_compute_thin_u(computationOptions);
424  m_computeFullV = ShouldComputeFullV || internal::should_svd_compute_full_v(computationOptions);
425  m_computeThinV = ShouldComputeThinV || internal::should_svd_compute_thin_v(computationOptions);
426 
427  eigen_assert(!(m_computeFullU && m_computeThinU) && "SVDBase: you can't ask for both full and thin U");
428  eigen_assert(!(m_computeFullV && m_computeThinV) && "SVDBase: you can't ask for both full and thin V");
429 
430  m_diagSize = (std::min)(m_rows, m_cols);
431  m_singularValues.resize(m_diagSize);
432  if(RowsAtCompileTime==Dynamic)
433  m_matrixU.resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0);
434  if(ColsAtCompileTime==Dynamic)
435  m_matrixV.resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0);
436 
437  return false;
438 }
439 
440 }// end namespace
441 
442 #endif // EIGEN_SVDBASE_H
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:182
Base class of SVD algorithms.
Definition: SVDBase.h:120
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SVDBase.h:310
Derived & setThreshold(const RealScalar &threshold)
Definition: SVDBase.h:247
const MatrixVType & matrixV() const
Definition: SVDBase.h:191
Derived & setThreshold(Default_t)
Definition: SVDBase.h:262
Index rank() const
Definition: SVDBase.h:222
const SingularValuesType & singularValues() const
Definition: SVDBase.h:203
bool computeV() const
Definition: SVDBase.h:284
Eigen::Index Index
Definition: SVDBase.h:130
bool computeU() const
Definition: SVDBase.h:282
RealScalar threshold() const
Definition: SVDBase.h:272
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
SVDBase()
Default Constructor.
Definition: SVDBase.h:358
const MatrixUType & matrixU() const
Definition: SVDBase.h:175
Index nonzeroSingularValues() const
Definition: SVDBase.h:210
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:71
ComputationInfo
Definition: Constants.h:442
@ NoQRPreconditioner
Definition: Constants.h:429
@ HouseholderQRPreconditioner
Definition: Constants.h:431
@ ColPivHouseholderQRPreconditioner
Definition: Constants.h:427
@ FullPivHouseholderQRPreconditioner
Definition: Constants.h:433
@ Success
Definition: Constants.h:444
@ ComputeFullV
Definition: Constants.h:399
@ ComputeThinV
Definition: Constants.h:401
@ ComputeFullU
Definition: Constants.h:395
@ ComputeThinU
Definition: Constants.h:397
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
const int Dynamic
Definition: Constants.h:24
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231