Medial Code Documentation
Loading...
Searching...
No Matches
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
19namespace Eigen {
47template<typename Derived>
49{
50
51public:
52 typedef typename internal::traits<Derived>::MatrixType MatrixType;
53 typedef typename MatrixType::Scalar Scalar;
54 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
55 typedef typename MatrixType::StorageIndex StorageIndex;
56 typedef Eigen::Index Index;
57 enum {
58 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
59 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
60 DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
61 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
62 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
63 MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
64 MatrixOptions = MatrixType::Options
65 };
66
69 typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
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
112 {
113 eigen_assert(m_isInitialized && "SVD is not initialized.");
114 return m_singularValues;
115 }
116
119 {
120 eigen_assert(m_isInitialized && "SVD is not initialized.");
121 return m_nonzeroSingularValues;
122 }
123
130 inline Index rank() const
131 {
132 using std::abs;
133 using std::max;
134 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
135 if(m_singularValues.size()==0) return 0;
136 RealScalar premultiplied_threshold = (max)(m_singularValues.coeff(0) * threshold(), (std::numeric_limits<RealScalar>::min)());
137 Index i = m_nonzeroSingularValues-1;
138 while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
139 return i+1;
140 }
141
157 {
158 m_usePrescribedThreshold = true;
159 m_prescribedThreshold = threshold;
160 return derived();
161 }
162
171 Derived& setThreshold(Default_t)
172 {
173 m_usePrescribedThreshold = false;
174 return derived();
175 }
176
182 {
183 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
184 return m_usePrescribedThreshold ? m_prescribedThreshold
185 : (std::max<Index>)(1,m_diagSize)*NumTraits<Scalar>::epsilon();
186 }
187
189 inline bool computeU() const { return m_computeFullU || m_computeThinU; }
191 inline bool computeV() const { return m_computeFullV || m_computeThinV; }
192
193 inline Index rows() const { return m_rows; }
194 inline Index cols() const { return m_cols; }
195
205 template<typename Rhs>
206 inline const Solve<Derived, Rhs>
207 solve(const MatrixBase<Rhs>& b) const
208 {
209 eigen_assert(m_isInitialized && "SVD is not initialized.");
210 eigen_assert(computeU() && computeV() && "SVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
211 return Solve<Derived, Rhs>(derived(), b.derived());
212 }
213
214 #ifndef EIGEN_PARSED_BY_DOXYGEN
215 template<typename RhsType, typename DstType>
217 void _solve_impl(const RhsType &rhs, DstType &dst) const;
218 #endif
219
220protected:
221
222 static void check_template_parameters()
223 {
224 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
225 }
226
227 // return true if already allocated
228 bool allocate(Index rows, Index cols, unsigned int computationOptions) ;
229
230 MatrixUType m_matrixU;
231 MatrixVType m_matrixV;
232 SingularValuesType m_singularValues;
233 bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
234 bool m_computeFullU, m_computeThinU;
235 bool m_computeFullV, m_computeThinV;
236 unsigned int m_computationOptions;
237 Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
238 RealScalar m_prescribedThreshold;
239
245 : m_isInitialized(false),
246 m_isAllocated(false),
247 m_usePrescribedThreshold(false),
248 m_computationOptions(0),
249 m_rows(-1), m_cols(-1), m_diagSize(0)
250 {
251 check_template_parameters();
252 }
253
254
255};
256
257#ifndef EIGEN_PARSED_BY_DOXYGEN
258template<typename Derived>
259template<typename RhsType, typename DstType>
260void SVDBase<Derived>::_solve_impl(const RhsType &rhs, DstType &dst) const
261{
262 eigen_assert(rhs.rows() == rows());
263
264 // A = U S V^*
265 // So A^{-1} = V S^{-1} U^*
266
268 Index l_rank = rank();
269 tmp.noalias() = m_matrixU.leftCols(l_rank).adjoint() * rhs;
270 tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
271 dst = m_matrixV.leftCols(l_rank) * tmp;
272}
273#endif
274
275template<typename MatrixType>
276bool SVDBase<MatrixType>::allocate(Index rows, Index cols, unsigned int computationOptions)
277{
278 eigen_assert(rows >= 0 && cols >= 0);
279
280 if (m_isAllocated &&
281 rows == m_rows &&
282 cols == m_cols &&
283 computationOptions == m_computationOptions)
284 {
285 return true;
286 }
287
288 m_rows = rows;
289 m_cols = cols;
290 m_isInitialized = false;
291 m_isAllocated = true;
292 m_computationOptions = computationOptions;
293 m_computeFullU = (computationOptions & ComputeFullU) != 0;
294 m_computeThinU = (computationOptions & ComputeThinU) != 0;
295 m_computeFullV = (computationOptions & ComputeFullV) != 0;
296 m_computeThinV = (computationOptions & ComputeThinV) != 0;
297 eigen_assert(!(m_computeFullU && m_computeThinU) && "SVDBase: you can't ask for both full and thin U");
298 eigen_assert(!(m_computeFullV && m_computeThinV) && "SVDBase: you can't ask for both full and thin V");
299 eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
300 "SVDBase: thin U and V are only available when your matrix has a dynamic number of columns.");
301
302 m_diagSize = (std::min)(m_rows, m_cols);
303 m_singularValues.resize(m_diagSize);
304 if(RowsAtCompileTime==Dynamic)
305 m_matrixU.resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0);
306 if(ColsAtCompileTime==Dynamic)
307 m_matrixV.resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0);
308
309 return false;
310}
311
312}// end namespace
313
314#endif // EIGEN_SVDBASE_H
Base class of SVD algorithms.
Definition SVDBase.h:49
Derived & setThreshold(const RealScalar &threshold)
Allows to prescribe a threshold to be used by certain methods, such as rank() and solve(),...
Definition SVDBase.h:156
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SVDBase.h:207
Index rank() const
Definition SVDBase.h:130
bool computeV() const
Definition SVDBase.h:191
Eigen::Index Index
Definition SVDBase.h:56
bool computeU() const
Definition SVDBase.h:189
Derived & setThreshold(Default_t)
Allows to come back to the default behavior, letting Eigen use its default formula for determining th...
Definition SVDBase.h:171
RealScalar threshold() const
Returns the threshold that will be used by certain methods such as rank().
Definition SVDBase.h:181
SVDBase()
Default Constructor.
Definition SVDBase.h:244
const SingularValuesType & singularValues() const
Definition SVDBase.h:111
const MatrixUType & matrixU() const
Definition SVDBase.h:83
const MatrixVType & matrixV() const
Definition SVDBase.h:99
Index nonzeroSingularValues() const
Definition SVDBase.h:118
Pseudo expression representing a solving operation.
Definition Solve.h:63
@ ComputeFullV
Used in JacobiSVD to indicate that the square matrix V is to be computed.
Definition Constants.h:387
@ ComputeThinV
Used in JacobiSVD to indicate that the thin matrix V is to be computed.
Definition Constants.h:389
@ ComputeFullU
Used in JacobiSVD to indicate that the square matrix U is to be computed.
Definition Constants.h:383
@ ComputeThinU
Used in JacobiSVD to indicate that the thin matrix U is to be computed.
Definition Constants.h:385
Holds information about the various numeric (i.e.
Definition NumTraits.h:108
Definition ForwardDeclarations.h:17
Definition Meta.h:30
Definition inference.c:32