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 {
47 template<typename Derived>
48 class SVDBase
49 {
50 
51 public:
53  typedef typename MatrixType::Scalar Scalar;
55  typedef typename MatrixType::StorageIndex StorageIndex;
56  typedef Eigen::Index Index;
57  enum {
58  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
59  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
61  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
62  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
64  MatrixOptions = MatrixType::Options
65  };
66 
70 
71  Derived& derived() { return *static_cast<Derived*>(this); }
72  const Derived& derived() const { return *static_cast<const Derived*>(this); }
73 
83  const MatrixUType& matrixU() const
84  {
85  eigen_assert(m_isInitialized && "SVD is not initialized.");
86  eigen_assert(computeU() && "This SVD decomposition didn't compute U. Did you ask for it?");
87  return m_matrixU;
88  }
89 
99  const MatrixVType& matrixV() const
100  {
101  eigen_assert(m_isInitialized && "SVD is not initialized.");
102  eigen_assert(computeV() && "This SVD decomposition didn't compute V. Did you ask for it?");
103  return m_matrixV;
104  }
105 
111  const SingularValuesType& singularValues() const
112  {
113  eigen_assert(m_isInitialized && "SVD is not initialized.");
114  return m_singularValues;
115  }
116 
118  Index nonzeroSingularValues() const
119  {
120  eigen_assert(m_isInitialized && "SVD is not initialized.");
122  }
123 
130  inline Index rank() const
131  {
132  using std::abs;
133  eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
134  if(m_singularValues.size()==0) return 0;
135  RealScalar premultiplied_threshold = numext::maxi<RealScalar>(m_singularValues.coeff(0) * threshold(), (std::numeric_limits<RealScalar>::min)());
136  Index i = m_nonzeroSingularValues-1;
137  while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
138  return i+1;
139  }
140 
155  Derived& setThreshold(const RealScalar& threshold)
156  {
159  return derived();
160  }
161 
171  {
172  m_usePrescribedThreshold = false;
173  return derived();
174  }
175 
180  RealScalar threshold() const
181  {
184  : (std::max<Index>)(1,m_diagSize)*NumTraits<Scalar>::epsilon();
185  }
186 
188  inline bool computeU() const { return m_computeFullU || m_computeThinU; }
190  inline bool computeV() const { return m_computeFullV || m_computeThinV; }
191 
192  inline Index rows() const { return m_rows; }
193  inline Index cols() const { return m_cols; }
194 
204  template<typename Rhs>
205  inline const Solve<Derived, Rhs>
206  solve(const MatrixBase<Rhs>& b) const
207  {
208  eigen_assert(m_isInitialized && "SVD is not initialized.");
209  eigen_assert(computeU() && computeV() && "SVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
210  return Solve<Derived, Rhs>(derived(), b.derived());
211  }
212 
213  #ifndef EIGEN_PARSED_BY_DOXYGEN
214  template<typename RhsType, typename DstType>
215  EIGEN_DEVICE_FUNC
216  void _solve_impl(const RhsType &rhs, DstType &dst) const;
217  #endif
218 
219 protected:
220 
222  {
224  }
225 
226  // return true if already allocated
227  bool allocate(Index rows, Index cols, unsigned int computationOptions) ;
228 
229  MatrixUType m_matrixU;
230  MatrixVType m_matrixV;
231  SingularValuesType m_singularValues;
235  unsigned int m_computationOptions;
238 
244  : m_isInitialized(false),
245  m_isAllocated(false),
246  m_usePrescribedThreshold(false),
247  m_computationOptions(0),
248  m_rows(-1), m_cols(-1), m_diagSize(0)
249  {
251  }
252 
253 
254 };
255 
256 #ifndef EIGEN_PARSED_BY_DOXYGEN
257 template<typename Derived>
258 template<typename RhsType, typename DstType>
259 void SVDBase<Derived>::_solve_impl(const RhsType &rhs, DstType &dst) const
260 {
261  eigen_assert(rhs.rows() == rows());
262 
263  // A = U S V^*
264  // So A^{-1} = V S^{-1} U^*
265 
267  Index l_rank = rank();
268  tmp.noalias() = m_matrixU.leftCols(l_rank).adjoint() * rhs;
269  tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
270  dst = m_matrixV.leftCols(l_rank) * tmp;
271 }
272 #endif
273 
274 template<typename MatrixType>
275 bool SVDBase<MatrixType>::allocate(Index rows, Index cols, unsigned int computationOptions)
276 {
277  eigen_assert(rows >= 0 && cols >= 0);
278 
279  if (m_isAllocated &&
280  rows == m_rows &&
281  cols == m_cols &&
282  computationOptions == m_computationOptions)
283  {
284  return true;
285  }
286 
287  m_rows = rows;
288  m_cols = cols;
289  m_isInitialized = false;
290  m_isAllocated = true;
291  m_computationOptions = computationOptions;
292  m_computeFullU = (computationOptions & ComputeFullU) != 0;
293  m_computeThinU = (computationOptions & ComputeThinU) != 0;
294  m_computeFullV = (computationOptions & ComputeFullV) != 0;
295  m_computeThinV = (computationOptions & ComputeThinV) != 0;
296  eigen_assert(!(m_computeFullU && m_computeThinU) && "SVDBase: you can't ask for both full and thin U");
297  eigen_assert(!(m_computeFullV && m_computeThinV) && "SVDBase: you can't ask for both full and thin V");
298  eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
299  "SVDBase: thin U and V are only available when your matrix has a dynamic number of columns.");
300 
301  m_diagSize = (std::min)(m_rows, m_cols);
307 
308  return false;
309 }
310 
311 }// end namespace
312 
313 #endif // EIGEN_SVDBASE_H
Derived & derived()
Definition: SVDBase.h:71
Derived & setThreshold(Default_t)
Definition: SVDBase.h:170
unsigned int m_computationOptions
Definition: SVDBase.h:235
Matrix< Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime > MatrixVType
Definition: SVDBase.h:68
const SingularValuesType & singularValues() const
Definition: SVDBase.h:111
bool m_computeFullV
Definition: SVDBase.h:234
Eigen::Index Index
Definition: SVDBase.h:56
MatrixType::StorageIndex StorageIndex
Definition: SVDBase.h:55
RealScalar threshold() const
Definition: SVDBase.h:180
MatrixUType m_matrixU
Definition: SVDBase.h:229
Definition: LDLT.h:16
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
Default_t
Definition: Constants.h:352
Index rank() const
Definition: SVDBase.h:130
NumTraits< typename MatrixType::Scalar >::Real RealScalar
Definition: SVDBase.h:54
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
EIGEN_DEVICE_FUNC void _solve_impl(const RhsType &rhs, DstType &dst) const
#define EIGEN_IMPLIES(a, b)
Definition: Macros.h:899
#define EIGEN_SIZE_MIN_PREFER_FIXED(a, b)
Definition: Macros.h:886
bool allocate(Index rows, Index cols, unsigned int computationOptions)
Definition: SVDBase.h:275
Index rows() const
Definition: SVDBase.h:192
bool computeV() const
Definition: SVDBase.h:190
bool m_isInitialized
Definition: SVDBase.h:232
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
bool m_usePrescribedThreshold
Definition: SVDBase.h:232
Base class of SVD algorithms.
Definition: SVDBase.h:48
MatrixVType m_matrixV
Definition: SVDBase.h:230
bool m_computeThinV
Definition: SVDBase.h:234
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Index m_cols
Definition: SVDBase.h:236
SVDBase()
Default Constructor.
Definition: SVDBase.h:243
#define eigen_assert(x)
Definition: Macros.h:577
Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime > MatrixUType
Definition: SVDBase.h:67
SingularValuesType m_singularValues
Definition: SVDBase.h:231
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:182
const MatrixVType & matrixV() const
Definition: SVDBase.h:99
Index m_rows
Definition: SVDBase.h:236
RealScalar m_prescribedThreshold
Definition: SVDBase.h:237
const MatrixUType & matrixU() const
Definition: SVDBase.h:83
Derived & setThreshold(const RealScalar &threshold)
Definition: SVDBase.h:155
const Derived & derived() const
Definition: SVDBase.h:72
internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
Definition: SVDBase.h:69
static void check_template_parameters()
Definition: SVDBase.h:221
bool computeU() const
Definition: SVDBase.h:188
Index nonzeroSingularValues() const
Definition: SVDBase.h:118
Index m_nonzeroSingularValues
Definition: SVDBase.h:236
bool m_computeFullU
Definition: SVDBase.h:233
Index m_diagSize
Definition: SVDBase.h:236
bool m_computeThinU
Definition: SVDBase.h:233
const int Dynamic
Definition: Constants.h:21
Pseudo expression representing a solving operation.
Definition: Solve.h:62
#define EIGEN_SIZE_MIN_PREFER_DYNAMIC(a, b)
Definition: Macros.h:878
EIGEN_DEVICE_FUNC const Scalar & b
internal::traits< Derived >::MatrixType MatrixType
Definition: SVDBase.h:52
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SVDBase.h:206
Index cols() const
Definition: SVDBase.h:193
bool m_isAllocated
Definition: SVDBase.h:232
MatrixType::Scalar Scalar
Definition: SVDBase.h:53


hebiros
Author(s): Xavier Artache , Matthew Tesch
autogenerated on Thu Sep 3 2020 04:09:07