16 #ifndef EIGEN_SVDBASE_H 17 #define EIGEN_SVDBASE_H 47 template<
typename Derived>
52 typedef typename internal::traits<Derived>::MatrixType MatrixType;
53 typedef typename MatrixType::Scalar Scalar;
55 typedef typename MatrixType::StorageIndex StorageIndex;
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
69 typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
71 Derived& derived() {
return *
static_cast<Derived*
>(
this); }
72 const Derived& derived()
const {
return *
static_cast<const Derived*
>(
this); }
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?");
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?");
113 eigen_assert(m_isInitialized &&
"SVD is not initialized.");
114 return m_singularValues;
120 eigen_assert(m_isInitialized &&
"SVD is not initialized.");
121 return m_nonzeroSingularValues;
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;
157 m_usePrescribedThreshold =
true;
172 m_usePrescribedThreshold =
false;
182 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
184 Index diagSize = (std::max<Index>)(1,m_diagSize);
185 return m_usePrescribedThreshold ? m_prescribedThreshold
190 inline bool computeU()
const {
return m_computeFullU || m_computeThinU; }
192 inline bool computeV()
const {
return m_computeFullV || m_computeThinV; }
194 inline Index rows()
const {
return m_rows; }
195 inline Index cols()
const {
return m_cols; }
206 template<
typename Rhs>
210 eigen_assert(m_isInitialized &&
"SVD is not initialized.");
211 eigen_assert(
computeU() &&
computeV() &&
"SVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
215 #ifndef EIGEN_PARSED_BY_DOXYGEN 216 template<
typename RhsType,
typename DstType>
218 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
223 static void check_template_parameters()
225 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
229 bool allocate(Index rows, Index cols,
unsigned int computationOptions) ;
231 MatrixUType m_matrixU;
232 MatrixVType m_matrixV;
233 SingularValuesType m_singularValues;
234 bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
235 bool m_computeFullU, m_computeThinU;
236 bool m_computeFullV, m_computeThinV;
237 unsigned int m_computationOptions;
238 Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
239 RealScalar m_prescribedThreshold;
246 : m_isInitialized(false),
247 m_isAllocated(false),
248 m_usePrescribedThreshold(false),
249 m_computationOptions(0),
250 m_rows(-1), m_cols(-1), m_diagSize(0)
252 check_template_parameters();
258 #ifndef EIGEN_PARSED_BY_DOXYGEN 259 template<
typename Derived>
260 template<
typename RhsType,
typename DstType>
263 eigen_assert(rhs.rows() == rows());
270 tmp.noalias() = m_matrixU.leftCols(l_rank).adjoint() * rhs;
271 tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
272 dst = m_matrixV.leftCols(l_rank) * tmp;
276 template<
typename MatrixType>
279 eigen_assert(rows >= 0 && cols >= 0);
284 computationOptions == m_computationOptions)
291 m_isInitialized =
false;
292 m_isAllocated =
true;
293 m_computationOptions = computationOptions;
294 m_computeFullU = (computationOptions &
ComputeFullU) != 0;
295 m_computeThinU = (computationOptions &
ComputeThinU) != 0;
296 m_computeFullV = (computationOptions &
ComputeFullV) != 0;
297 m_computeThinV = (computationOptions &
ComputeThinV) != 0;
298 eigen_assert(!(m_computeFullU && m_computeThinU) &&
"SVDBase: you can't ask for both full and thin U");
299 eigen_assert(!(m_computeFullV && m_computeThinV) &&
"SVDBase: you can't ask for both full and thin V");
300 eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==
Dynamic) &&
301 "SVDBase: thin U and V are only available when your matrix has a dynamic number of columns.");
303 m_diagSize = (std::min)(m_rows, m_cols);
304 m_singularValues.resize(m_diagSize);
306 m_matrixU.
resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0);
308 m_matrixV.
resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0);
315 #endif // EIGEN_SVDBASE_H Definition: Constants.h:383
Index rank() const
Definition: SVDBase.h:130
Derived & setThreshold(Default_t)
Definition: SVDBase.h:170
Definition: Constants.h:389
const MatrixUType & matrixU() const
Definition: SVDBase.h:83
Eigen::Index Index
Definition: SVDBase.h:56
Namespace containing all symbols from the Eigen library.
Definition: Core:306
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
Derived & derived()
Definition: EigenBase.h:45
void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:279
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SVDBase.h:208
bool computeV() const
Definition: SVDBase.h:192
Base class of SVD algorithms.
Definition: SVDBase.h:48
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
SVDBase()
Default Constructor.
Definition: SVDBase.h:245
Derived & setThreshold(const RealScalar &threshold)
Definition: SVDBase.h:155
RealScalar threshold() const
Definition: SVDBase.h:180
Index nonzeroSingularValues() const
Definition: SVDBase.h:118
const SingularValuesType & singularValues() const
Definition: SVDBase.h:111
bool computeU() const
Definition: SVDBase.h:190
const MatrixVType & matrixV() const
Definition: SVDBase.h:99
const int Dynamic
Definition: Constants.h:21
Pseudo expression representing a solving operation.
Definition: Solve.h:62
Definition: Constants.h:387
Definition: Constants.h:385
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48