Eigen  3.3.0
 
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;
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
111 const SingularValuesType& singularValues() const
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 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 {
157 m_usePrescribedThreshold = true;
158 m_prescribedThreshold = threshold;
159 return derived();
160 }
161
170 Derived& setThreshold(Default_t)
171 {
172 m_usePrescribedThreshold = false;
173 return derived();
174 }
175
180 RealScalar threshold() const
181 {
182 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
183 return m_usePrescribedThreshold ? m_prescribedThreshold
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
219protected:
220
221 static void check_template_parameters()
222 {
223 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
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;
232 bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
233 bool m_computeFullU, m_computeThinU;
234 bool m_computeFullV, m_computeThinV;
235 unsigned int m_computationOptions;
236 Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
237 RealScalar m_prescribedThreshold;
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 {
250 check_template_parameters();
251 }
252
253
254};
255
256#ifndef EIGEN_PARSED_BY_DOXYGEN
257template<typename Derived>
258template<typename RhsType, typename DstType>
259void 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
266 Matrix<Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp;
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
274template<typename MatrixType>
275bool 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);
302 m_singularValues.resize(m_diagSize);
303 if(RowsAtCompileTime==Dynamic)
304 m_matrixU.resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0);
305 if(ColsAtCompileTime==Dynamic)
306 m_matrixV.resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0);
307
308 return false;
309}
310
311}// end namespace
312
313#endif // EIGEN_SVDBASE_H
Derived & derived()
Definition: EigenBase.h:44
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:180
Base class of SVD algorithms.
Definition: SVDBase.h:49
Derived & setThreshold(const RealScalar &threshold)
Definition: SVDBase.h:155
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SVDBase.h:206
Index rank() const
Definition: SVDBase.h:130
bool computeV() const
Definition: SVDBase.h:190
Eigen::Index Index
Definition: SVDBase.h:56
bool computeU() const
Definition: SVDBase.h:188
Derived & setThreshold(Default_t)
Definition: SVDBase.h:170
RealScalar threshold() const
Definition: SVDBase.h:180
SVDBase()
Default Constructor.
Definition: SVDBase.h:243
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
Definition: Constants.h:387
@ ComputeThinV
Definition: Constants.h:389
@ ComputeFullU
Definition: Constants.h:383
@ ComputeThinU
Definition: Constants.h:385
Namespace containing all symbols from the Eigen library.
Definition: Core:287
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:151