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 namespace Eigen {
20 
21 namespace internal {
22 template<typename Derived> struct traits<SVDBase<Derived> >
23  : traits<Derived>
24 {
25  typedef MatrixXpr XprKind;
27  typedef int StorageIndex;
28  enum { Flags = 0 };
29 };
30 }
31 
62 template<typename Derived> class SVDBase
63  : public SolverBase<SVDBase<Derived> >
64 {
65 public:
66 
67  template<typename Derived_>
69 
71  typedef typename MatrixType::Scalar Scalar;
74  typedef Eigen::Index Index;
75  enum {
76  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
77  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
79  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
80  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
82  MatrixOptions = MatrixType::Options
83  };
84 
88 
89  Derived& derived() { return *static_cast<Derived*>(this); }
90  const Derived& derived() const { return *static_cast<const Derived*>(this); }
91 
101  const MatrixUType& matrixU() const
102  {
104  eigen_assert(computeU() && "This SVD decomposition didn't compute U. Did you ask for it?");
105  return m_matrixU;
106  }
107 
117  const MatrixVType& matrixV() const
118  {
120  eigen_assert(computeV() && "This SVD decomposition didn't compute V. Did you ask for it?");
121  return m_matrixV;
122  }
123 
130  {
132  return m_singularValues;
133  }
134 
137  {
140  }
141 
148  inline Index rank() const
149  {
150  using std::abs;
152  if(m_singularValues.size()==0) return 0;
153  RealScalar premultiplied_threshold = numext::maxi<RealScalar>(m_singularValues.coeff(0) * threshold(), (std::numeric_limits<RealScalar>::min)());
155  while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
156  return i+1;
157  }
158 
174  {
177  return derived();
178  }
179 
189  {
190  m_usePrescribedThreshold = false;
191  return derived();
192  }
193 
199  {
201  // this temporary is needed to workaround a MSVC issue
202  Index diagSize = (std::max<Index>)(1,m_diagSize);
205  }
206 
208  inline bool computeU() const { return m_computeFullU || m_computeThinU; }
210  inline bool computeV() const { return m_computeFullV || m_computeThinV; }
211 
212  inline Index rows() const { return m_rows; }
213  inline Index cols() const { return m_cols; }
214 
215  #ifdef EIGEN_PARSED_BY_DOXYGEN
216 
225  template<typename Rhs>
226  inline const Solve<Derived, Rhs>
227  solve(const MatrixBase<Rhs>& b) const;
228  #endif
229 
230 
237  {
238  eigen_assert(m_isInitialized && "SVD is not initialized.");
239  return m_info;
240  }
241 
242  #ifndef EIGEN_PARSED_BY_DOXYGEN
243  template<typename RhsType, typename DstType>
244  void _solve_impl(const RhsType &rhs, DstType &dst) const;
245 
246  template<bool Conjugate, typename RhsType, typename DstType>
247  void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
248  #endif
249 
250 protected:
251 
253  {
255  }
256 
258  eigen_assert(m_isInitialized && "SVD is not initialized.");
259  }
260 
261  template<bool Transpose_, typename Rhs>
262  void _check_solve_assertion(const Rhs& b) const {
265  eigen_assert(computeU() && computeV() && "SVDBase::solve(): Both unitaries U and V are required to be computed (thin unitaries suffice).");
266  eigen_assert((Transpose_?cols():rows())==b.rows() && "SVDBase::solve(): invalid number of rows of the right hand side matrix b");
267  }
268 
269  // return true if already allocated
270  bool allocate(Index rows, Index cols, unsigned int computationOptions) ;
271 
279  unsigned int m_computationOptions;
282 
288  : m_info(Success),
289  m_isInitialized(false),
290  m_isAllocated(false),
292  m_computeFullU(false),
293  m_computeThinU(false),
294  m_computeFullV(false),
295  m_computeThinV(false),
297  m_rows(-1), m_cols(-1), m_diagSize(0)
298  {
300  }
301 
302 
303 };
304 
305 #ifndef EIGEN_PARSED_BY_DOXYGEN
306 template<typename Derived>
307 template<typename RhsType, typename DstType>
308 void SVDBase<Derived>::_solve_impl(const RhsType &rhs, DstType &dst) const
309 {
310  // A = U S V^*
311  // So A^{-1} = V S^{-1} U^*
312 
314  Index l_rank = rank();
315  tmp.noalias() = m_matrixU.leftCols(l_rank).adjoint() * rhs;
316  tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
317  dst = m_matrixV.leftCols(l_rank) * tmp;
318 }
319 
320 template<typename Derived>
321 template<bool Conjugate, typename RhsType, typename DstType>
322 void SVDBase<Derived>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
323 {
324  // A = U S V^*
325  // So A^{-*} = U S^{-1} V^*
326  // And A^{-T} = U_conj S^{-1} V^T
328  Index l_rank = rank();
329 
330  tmp.noalias() = m_matrixV.leftCols(l_rank).transpose().template conjugateIf<Conjugate>() * rhs;
331  tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
332  dst = m_matrixU.template conjugateIf<!Conjugate>().leftCols(l_rank) * tmp;
333 }
334 #endif
335 
336 template<typename MatrixType>
337 bool SVDBase<MatrixType>::allocate(Index rows, Index cols, unsigned int computationOptions)
338 {
339  eigen_assert(rows >= 0 && cols >= 0);
340 
341  if (m_isAllocated &&
342  rows == m_rows &&
343  cols == m_cols &&
344  computationOptions == m_computationOptions)
345  {
346  return true;
347  }
348 
349  m_rows = rows;
350  m_cols = cols;
351  m_info = Success;
352  m_isInitialized = false;
353  m_isAllocated = true;
354  m_computationOptions = computationOptions;
355  m_computeFullU = (computationOptions & ComputeFullU) != 0;
356  m_computeThinU = (computationOptions & ComputeThinU) != 0;
357  m_computeFullV = (computationOptions & ComputeFullV) != 0;
358  m_computeThinV = (computationOptions & ComputeThinV) != 0;
359  eigen_assert(!(m_computeFullU && m_computeThinU) && "SVDBase: you can't ask for both full and thin U");
360  eigen_assert(!(m_computeFullV && m_computeThinV) && "SVDBase: you can't ask for both full and thin V");
361  eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
362  "SVDBase: thin U and V are only available when your matrix has a dynamic number of columns.");
363 
364  m_diagSize = (std::min)(m_rows, m_cols);
365  m_singularValues.resize(m_diagSize);
366  if(RowsAtCompileTime==Dynamic)
367  m_matrixU.resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0);
368  if(ColsAtCompileTime==Dynamic)
369  m_matrixV.resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0);
370 
371  return false;
372 }
373 
374 }// end namespace
375 
376 #endif // EIGEN_SVDBASE_H
Eigen::SVDBase::SVDBase
SVDBase()
Default Constructor.
Definition: SVDBase.h:287
Eigen::SVDBase::m_computeFullU
bool m_computeFullU
Definition: SVDBase.h:277
Eigen::SVDBase::m_isAllocated
bool m_isAllocated
Definition: SVDBase.h:276
Eigen::SVDBase::_check_solve_assertion
void _check_solve_assertion(const Rhs &b) const
Definition: SVDBase.h:262
Eigen::SVDBase::m_cols
Index m_cols
Definition: SVDBase.h:280
Eigen::MatrixXpr
Definition: Constants.h:522
EIGEN_DEVICE_FUNC
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:976
Eigen
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
Eigen::SVDBase::m_usePrescribedThreshold
bool m_usePrescribedThreshold
Definition: SVDBase.h:276
Eigen::SVDBase::m_computeThinU
bool m_computeThinU
Definition: SVDBase.h:277
Eigen::ComputeFullV
@ ComputeFullV
Definition: Constants.h:397
Eigen::internal::traits< SVDBase< Derived > >::StorageKind
SolverStorage StorageKind
Definition: SVDBase.h:26
Eigen::SVDBase
Base class of SVD algorithms.
Definition: ForwardDeclarations.h:277
b
Scalar * b
Definition: benchVecAdd.cpp:17
eigen_assert
#define eigen_assert(x)
Definition: Macros.h:1037
Eigen::SVDBase::allocate
bool allocate(Index rows, Index cols, unsigned int computationOptions)
Definition: SVDBase.h:337
Eigen::SVDBase::setThreshold
Derived & setThreshold(const RealScalar &threshold)
Definition: SVDBase.h:173
Eigen::SVDBase::threshold
RealScalar threshold() const
Definition: SVDBase.h:198
Eigen::ComputeFullU
@ ComputeFullU
Definition: Constants.h:393
Eigen::Success
@ Success
Definition: Constants.h:442
Eigen::SVDBase::MatrixVType
Matrix< Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime > MatrixVType
Definition: SVDBase.h:86
Eigen::SVDBase::MatrixType
internal::traits< Derived >::MatrixType MatrixType
Definition: SVDBase.h:70
Eigen::SVDBase::m_singularValues
SingularValuesType m_singularValues
Definition: SVDBase.h:274
Eigen::internal::solve_assertion
Definition: SolverBase.h:18
Eigen::SVDBase::m_computeFullV
bool m_computeFullV
Definition: SVDBase.h:278
Eigen::SVDBase::derived
const Derived & derived() const
Definition: SVDBase.h:90
Eigen::SVDBase::MatrixOptions
@ MatrixOptions
Definition: SVDBase.h:82
Eigen::SVDBase::m_computeThinV
bool m_computeThinV
Definition: SVDBase.h:278
rows
int rows
Definition: Tutorial_commainit_02.cpp:1
Eigen::SVDBase::RowsAtCompileTime
@ RowsAtCompileTime
Definition: SVDBase.h:76
EIGEN_ONLY_USED_FOR_DEBUG
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
Definition: Macros.h:1049
Eigen::SVDBase::StorageIndex
Eigen::internal::traits< SVDBase >::StorageIndex StorageIndex
Definition: SVDBase.h:73
Eigen::SVDBase::RealScalar
NumTraits< typename MatrixType::Scalar >::Real RealScalar
Definition: SVDBase.h:72
Eigen::internal::traits< SVDBase< Derived > >::StorageIndex
int StorageIndex
Definition: SVDBase.h:27
epsilon
static double epsilon
Definition: testRot3.cpp:37
Eigen::internal::true_type
Definition: Meta.h:96
Eigen::SVDBase::m_info
ComputationInfo m_info
Definition: SVDBase.h:275
Eigen::SVDBase::_solve_impl
void _solve_impl(const RhsType &rhs, DstType &dst) const
Definition: SVDBase.h:308
Eigen::SolverStorage
Definition: Constants.h:513
Eigen::ComputeThinU
@ ComputeThinU
Definition: Constants.h:395
Eigen::SVDBase::m_matrixU
MatrixUType m_matrixU
Definition: SVDBase.h:272
Eigen::Dynamic
const int Dynamic
Definition: Constants.h:22
Eigen::SVDBase::derived
Derived & derived()
Definition: SVDBase.h:89
Eigen::SVDBase::computeV
bool computeV() const
Definition: SVDBase.h:210
Eigen::SVDBase::rank
Index rank() const
Definition: SVDBase.h:148
EIGEN_SIZE_MIN_PREFER_DYNAMIC
#define EIGEN_SIZE_MIN_PREFER_DYNAMIC(a, b)
Definition: Macros.h:1294
Eigen::SVDBase::m_nonzeroSingularValues
Index m_nonzeroSingularValues
Definition: SVDBase.h:280
Eigen::SVDBase::nonzeroSingularValues
Index nonzeroSingularValues() const
Definition: SVDBase.h:136
Eigen::SVDBase::rows
Index rows() const
Definition: SVDBase.h:212
Eigen::SVDBase::cols
Index cols() const
Definition: SVDBase.h:213
Eigen::SVDBase::matrixV
const MatrixVType & matrixV() const
Definition: SVDBase.h:117
Eigen::SVDBase::computeU
bool computeU() const
Definition: SVDBase.h:208
Eigen::SVDBase::singularValues
const SingularValuesType & singularValues() const
Definition: SVDBase.h:129
Eigen::SVDBase::check_template_parameters
static void check_template_parameters()
Definition: SVDBase.h:252
Eigen::SVDBase::MaxDiagSizeAtCompileTime
@ MaxDiagSizeAtCompileTime
Definition: SVDBase.h:81
Eigen::SVDBase::MatrixUType
Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime > MatrixUType
Definition: SVDBase.h:85
Eigen::ComputeThinV
@ ComputeThinV
Definition: Constants.h:399
Eigen::Solve
Pseudo expression representing a solving operation.
Definition: Solve.h:62
Eigen::SVDBase::Index
Eigen::Index Index
Definition: SVDBase.h:74
Eigen::internal::traits
Definition: ForwardDeclarations.h:17
Eigen::SVDBase::SingularValuesType
internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
Definition: SVDBase.h:87
Eigen::SVDBase::ColsAtCompileTime
@ ColsAtCompileTime
Definition: SVDBase.h:77
Eigen::Default_t
Default_t
Definition: Constants.h:362
Eigen::internal::Rhs
@ Rhs
Definition: TensorContractionMapper.h:18
Eigen::SVDBase::MaxColsAtCompileTime
@ MaxColsAtCompileTime
Definition: SVDBase.h:80
Eigen::SVDBase::DiagSizeAtCompileTime
@ DiagSizeAtCompileTime
Definition: SVDBase.h:78
Eigen::internal::traits< SVDBase< Derived > >::XprKind
MatrixXpr XprKind
Definition: SVDBase.h:25
Eigen::SVDBase::m_diagSize
Index m_diagSize
Definition: SVDBase.h:280
Eigen::SVDBase::m_rows
Index m_rows
Definition: SVDBase.h:280
Eigen::SVDBase::_check_compute_assertions
void _check_compute_assertions() const
Definition: SVDBase.h:257
min
#define min(a, b)
Definition: datatypes.h:19
Eigen::Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime >
abs
#define abs(x)
Definition: datatypes.h:17
internal
Definition: BandTriangularSolver.h:13
Eigen::MatrixBase
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
EIGEN_IMPLIES
#define EIGEN_IMPLIES(a, b)
Definition: Macros.h:1315
Eigen::SVDBase::Scalar
MatrixType::Scalar Scalar
Definition: SVDBase.h:71
Eigen::SVDBase::m_isInitialized
bool m_isInitialized
Definition: SVDBase.h:276
Eigen::SVDBase::MaxRowsAtCompileTime
@ MaxRowsAtCompileTime
Definition: SVDBase.h:79
cols
int cols
Definition: Tutorial_commainit_02.cpp:1
Eigen::SVDBase::m_computationOptions
unsigned int m_computationOptions
Definition: SVDBase.h:279
Eigen::ComputationInfo
ComputationInfo
Definition: Constants.h:440
Eigen::SVDBase::info
EIGEN_DEVICE_FUNC ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SVDBase.h:236
Eigen::SVDBase::m_matrixV
MatrixVType m_matrixV
Definition: SVDBase.h:273
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:232
Eigen::SVDBase::_solve_impl_transposed
void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const
Definition: SVDBase.h:322
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
EIGEN_SIZE_MIN_PREFER_FIXED
#define EIGEN_SIZE_MIN_PREFER_FIXED(a, b)
Definition: Macros.h:1302
Eigen::SVDBase::setThreshold
Derived & setThreshold(Default_t)
Definition: SVDBase.h:188
Eigen::SVDBase::m_prescribedThreshold
RealScalar m_prescribedThreshold
Definition: SVDBase.h:281
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46
Eigen::SVDBase::matrixU
const MatrixUType & matrixU() const
Definition: SVDBase.h:101
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
EIGEN_STATIC_ASSERT_NON_INTEGER
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:187


gtsam
Author(s):
autogenerated on Thu Jun 13 2024 03:06:53