1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library 2*bf2c3715SXin Li // for linear algebra. 3*bf2c3715SXin Li // 4*bf2c3715SXin Li // Copyright (C) 2009-2010 Benoit Jacob <[email protected]> 5*bf2c3715SXin Li // Copyright (C) 2014 Gael Guennebaud <[email protected]> 6*bf2c3715SXin Li // 7*bf2c3715SXin Li // Copyright (C) 2013 Gauthier Brun <[email protected]> 8*bf2c3715SXin Li // Copyright (C) 2013 Nicolas Carre <[email protected]> 9*bf2c3715SXin Li // Copyright (C) 2013 Jean Ceccato <[email protected]> 10*bf2c3715SXin Li // Copyright (C) 2013 Pierre Zoppitelli <[email protected]> 11*bf2c3715SXin Li // 12*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla 13*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed 14*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 15*bf2c3715SXin Li 16*bf2c3715SXin Li #ifndef EIGEN_SVDBASE_H 17*bf2c3715SXin Li #define EIGEN_SVDBASE_H 18*bf2c3715SXin Li 19*bf2c3715SXin Li namespace Eigen { 20*bf2c3715SXin Li 21*bf2c3715SXin Li namespace internal { 22*bf2c3715SXin Li template<typename Derived> struct traits<SVDBase<Derived> > 23*bf2c3715SXin Li : traits<Derived> 24*bf2c3715SXin Li { 25*bf2c3715SXin Li typedef MatrixXpr XprKind; 26*bf2c3715SXin Li typedef SolverStorage StorageKind; 27*bf2c3715SXin Li typedef int StorageIndex; 28*bf2c3715SXin Li enum { Flags = 0 }; 29*bf2c3715SXin Li }; 30*bf2c3715SXin Li } 31*bf2c3715SXin Li 32*bf2c3715SXin Li /** \ingroup SVD_Module 33*bf2c3715SXin Li * 34*bf2c3715SXin Li * 35*bf2c3715SXin Li * \class SVDBase 36*bf2c3715SXin Li * 37*bf2c3715SXin Li * \brief Base class of SVD algorithms 38*bf2c3715SXin Li * 39*bf2c3715SXin Li * \tparam Derived the type of the actual SVD decomposition 40*bf2c3715SXin Li * 41*bf2c3715SXin Li * SVD decomposition consists in decomposing any n-by-p matrix \a A as a product 42*bf2c3715SXin Li * \f[ A = U S V^* \f] 43*bf2c3715SXin Li * where \a U is a n-by-n unitary, \a V is a p-by-p unitary, and \a S is a n-by-p real positive matrix which is zero outside of its main diagonal; 44*bf2c3715SXin Li * the diagonal entries of S are known as the \em singular \em values of \a A and the columns of \a U and \a V are known as the left 45*bf2c3715SXin Li * and right \em singular \em vectors of \a A respectively. 46*bf2c3715SXin Li * 47*bf2c3715SXin Li * Singular values are always sorted in decreasing order. 48*bf2c3715SXin Li * 49*bf2c3715SXin Li * 50*bf2c3715SXin Li * You can ask for only \em thin \a U or \a V to be computed, meaning the following. In case of a rectangular n-by-p matrix, letting \a m be the 51*bf2c3715SXin Li * smaller value among \a n and \a p, there are only \a m singular vectors; the remaining columns of \a U and \a V do not correspond to actual 52*bf2c3715SXin Li * singular vectors. Asking for \em thin \a U or \a V means asking for only their \a m first columns to be formed. So \a U is then a n-by-m matrix, 53*bf2c3715SXin Li * and \a V is then a p-by-m matrix. Notice that thin \a U and \a V are all you need for (least squares) solving. 54*bf2c3715SXin Li * 55*bf2c3715SXin Li * The status of the computation can be retrived using the \a info() method. Unless \a info() returns \a Success, the results should be not 56*bf2c3715SXin Li * considered well defined. 57*bf2c3715SXin Li * 58*bf2c3715SXin Li * If the input matrix has inf or nan coefficients, the result of the computation is undefined, and \a info() will return \a InvalidInput, but the computation is guaranteed to 59*bf2c3715SXin Li * terminate in finite (and reasonable) time. 60*bf2c3715SXin Li * \sa class BDCSVD, class JacobiSVD 61*bf2c3715SXin Li */ 62*bf2c3715SXin Li template<typename Derived> class SVDBase 63*bf2c3715SXin Li : public SolverBase<SVDBase<Derived> > 64*bf2c3715SXin Li { 65*bf2c3715SXin Li public: 66*bf2c3715SXin Li 67*bf2c3715SXin Li template<typename Derived_> 68*bf2c3715SXin Li friend struct internal::solve_assertion; 69*bf2c3715SXin Li 70*bf2c3715SXin Li typedef typename internal::traits<Derived>::MatrixType MatrixType; 71*bf2c3715SXin Li typedef typename MatrixType::Scalar Scalar; 72*bf2c3715SXin Li typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar; 73*bf2c3715SXin Li typedef typename Eigen::internal::traits<SVDBase>::StorageIndex StorageIndex; 74*bf2c3715SXin Li typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3 75*bf2c3715SXin Li enum { 76*bf2c3715SXin Li RowsAtCompileTime = MatrixType::RowsAtCompileTime, 77*bf2c3715SXin Li ColsAtCompileTime = MatrixType::ColsAtCompileTime, 78*bf2c3715SXin Li DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime), 79*bf2c3715SXin Li MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, 80*bf2c3715SXin Li MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime, 81*bf2c3715SXin Li MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime), 82*bf2c3715SXin Li MatrixOptions = MatrixType::Options 83*bf2c3715SXin Li }; 84*bf2c3715SXin Li 85*bf2c3715SXin Li typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixUType; 86*bf2c3715SXin Li typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime> MatrixVType; 87*bf2c3715SXin Li typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType; 88*bf2c3715SXin Li 89*bf2c3715SXin Li Derived& derived() { return *static_cast<Derived*>(this); } 90*bf2c3715SXin Li const Derived& derived() const { return *static_cast<const Derived*>(this); } 91*bf2c3715SXin Li 92*bf2c3715SXin Li /** \returns the \a U matrix. 93*bf2c3715SXin Li * 94*bf2c3715SXin Li * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p, 95*bf2c3715SXin Li * the U matrix is n-by-n if you asked for \link Eigen::ComputeFullU ComputeFullU \endlink, and is n-by-m if you asked for \link Eigen::ComputeThinU ComputeThinU \endlink. 96*bf2c3715SXin Li * 97*bf2c3715SXin Li * The \a m first columns of \a U are the left singular vectors of the matrix being decomposed. 98*bf2c3715SXin Li * 99*bf2c3715SXin Li * This method asserts that you asked for \a U to be computed. 100*bf2c3715SXin Li */ 101*bf2c3715SXin Li const MatrixUType& matrixU() const 102*bf2c3715SXin Li { 103*bf2c3715SXin Li _check_compute_assertions(); 104*bf2c3715SXin Li eigen_assert(computeU() && "This SVD decomposition didn't compute U. Did you ask for it?"); 105*bf2c3715SXin Li return m_matrixU; 106*bf2c3715SXin Li } 107*bf2c3715SXin Li 108*bf2c3715SXin Li /** \returns the \a V matrix. 109*bf2c3715SXin Li * 110*bf2c3715SXin Li * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p, 111*bf2c3715SXin Li * the V matrix is p-by-p if you asked for \link Eigen::ComputeFullV ComputeFullV \endlink, and is p-by-m if you asked for \link Eigen::ComputeThinV ComputeThinV \endlink. 112*bf2c3715SXin Li * 113*bf2c3715SXin Li * The \a m first columns of \a V are the right singular vectors of the matrix being decomposed. 114*bf2c3715SXin Li * 115*bf2c3715SXin Li * This method asserts that you asked for \a V to be computed. 116*bf2c3715SXin Li */ 117*bf2c3715SXin Li const MatrixVType& matrixV() const 118*bf2c3715SXin Li { 119*bf2c3715SXin Li _check_compute_assertions(); 120*bf2c3715SXin Li eigen_assert(computeV() && "This SVD decomposition didn't compute V. Did you ask for it?"); 121*bf2c3715SXin Li return m_matrixV; 122*bf2c3715SXin Li } 123*bf2c3715SXin Li 124*bf2c3715SXin Li /** \returns the vector of singular values. 125*bf2c3715SXin Li * 126*bf2c3715SXin Li * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p, the 127*bf2c3715SXin Li * returned vector has size \a m. Singular values are always sorted in decreasing order. 128*bf2c3715SXin Li */ 129*bf2c3715SXin Li const SingularValuesType& singularValues() const 130*bf2c3715SXin Li { 131*bf2c3715SXin Li _check_compute_assertions(); 132*bf2c3715SXin Li return m_singularValues; 133*bf2c3715SXin Li } 134*bf2c3715SXin Li 135*bf2c3715SXin Li /** \returns the number of singular values that are not exactly 0 */ 136*bf2c3715SXin Li Index nonzeroSingularValues() const 137*bf2c3715SXin Li { 138*bf2c3715SXin Li _check_compute_assertions(); 139*bf2c3715SXin Li return m_nonzeroSingularValues; 140*bf2c3715SXin Li } 141*bf2c3715SXin Li 142*bf2c3715SXin Li /** \returns the rank of the matrix of which \c *this is the SVD. 143*bf2c3715SXin Li * 144*bf2c3715SXin Li * \note This method has to determine which singular values should be considered nonzero. 145*bf2c3715SXin Li * For that, it uses the threshold value that you can control by calling 146*bf2c3715SXin Li * setThreshold(const RealScalar&). 147*bf2c3715SXin Li */ 148*bf2c3715SXin Li inline Index rank() const 149*bf2c3715SXin Li { 150*bf2c3715SXin Li using std::abs; 151*bf2c3715SXin Li _check_compute_assertions(); 152*bf2c3715SXin Li if(m_singularValues.size()==0) return 0; 153*bf2c3715SXin Li RealScalar premultiplied_threshold = numext::maxi<RealScalar>(m_singularValues.coeff(0) * threshold(), (std::numeric_limits<RealScalar>::min)()); 154*bf2c3715SXin Li Index i = m_nonzeroSingularValues-1; 155*bf2c3715SXin Li while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i; 156*bf2c3715SXin Li return i+1; 157*bf2c3715SXin Li } 158*bf2c3715SXin Li 159*bf2c3715SXin Li /** Allows to prescribe a threshold to be used by certain methods, such as rank() and solve(), 160*bf2c3715SXin Li * which need to determine when singular values are to be considered nonzero. 161*bf2c3715SXin Li * This is not used for the SVD decomposition itself. 162*bf2c3715SXin Li * 163*bf2c3715SXin Li * When it needs to get the threshold value, Eigen calls threshold(). 164*bf2c3715SXin Li * The default is \c NumTraits<Scalar>::epsilon() 165*bf2c3715SXin Li * 166*bf2c3715SXin Li * \param threshold The new value to use as the threshold. 167*bf2c3715SXin Li * 168*bf2c3715SXin Li * A singular value will be considered nonzero if its value is strictly greater than 169*bf2c3715SXin Li * \f$ \vert singular value \vert \leqslant threshold \times \vert max singular value \vert \f$. 170*bf2c3715SXin Li * 171*bf2c3715SXin Li * If you want to come back to the default behavior, call setThreshold(Default_t) 172*bf2c3715SXin Li */ 173*bf2c3715SXin Li Derived& setThreshold(const RealScalar& threshold) 174*bf2c3715SXin Li { 175*bf2c3715SXin Li m_usePrescribedThreshold = true; 176*bf2c3715SXin Li m_prescribedThreshold = threshold; 177*bf2c3715SXin Li return derived(); 178*bf2c3715SXin Li } 179*bf2c3715SXin Li 180*bf2c3715SXin Li /** Allows to come back to the default behavior, letting Eigen use its default formula for 181*bf2c3715SXin Li * determining the threshold. 182*bf2c3715SXin Li * 183*bf2c3715SXin Li * You should pass the special object Eigen::Default as parameter here. 184*bf2c3715SXin Li * \code svd.setThreshold(Eigen::Default); \endcode 185*bf2c3715SXin Li * 186*bf2c3715SXin Li * See the documentation of setThreshold(const RealScalar&). 187*bf2c3715SXin Li */ 188*bf2c3715SXin Li Derived& setThreshold(Default_t) 189*bf2c3715SXin Li { 190*bf2c3715SXin Li m_usePrescribedThreshold = false; 191*bf2c3715SXin Li return derived(); 192*bf2c3715SXin Li } 193*bf2c3715SXin Li 194*bf2c3715SXin Li /** Returns the threshold that will be used by certain methods such as rank(). 195*bf2c3715SXin Li * 196*bf2c3715SXin Li * See the documentation of setThreshold(const RealScalar&). 197*bf2c3715SXin Li */ 198*bf2c3715SXin Li RealScalar threshold() const 199*bf2c3715SXin Li { 200*bf2c3715SXin Li eigen_assert(m_isInitialized || m_usePrescribedThreshold); 201*bf2c3715SXin Li // this temporary is needed to workaround a MSVC issue 202*bf2c3715SXin Li Index diagSize = (std::max<Index>)(1,m_diagSize); 203*bf2c3715SXin Li return m_usePrescribedThreshold ? m_prescribedThreshold 204*bf2c3715SXin Li : RealScalar(diagSize)*NumTraits<Scalar>::epsilon(); 205*bf2c3715SXin Li } 206*bf2c3715SXin Li 207*bf2c3715SXin Li /** \returns true if \a U (full or thin) is asked for in this SVD decomposition */ 208*bf2c3715SXin Li inline bool computeU() const { return m_computeFullU || m_computeThinU; } 209*bf2c3715SXin Li /** \returns true if \a V (full or thin) is asked for in this SVD decomposition */ 210*bf2c3715SXin Li inline bool computeV() const { return m_computeFullV || m_computeThinV; } 211*bf2c3715SXin Li 212*bf2c3715SXin Li inline Index rows() const { return m_rows; } 213*bf2c3715SXin Li inline Index cols() const { return m_cols; } 214*bf2c3715SXin Li 215*bf2c3715SXin Li #ifdef EIGEN_PARSED_BY_DOXYGEN 216*bf2c3715SXin Li /** \returns a (least squares) solution of \f$ A x = b \f$ using the current SVD decomposition of A. 217*bf2c3715SXin Li * 218*bf2c3715SXin Li * \param b the right-hand-side of the equation to solve. 219*bf2c3715SXin Li * 220*bf2c3715SXin Li * \note Solving requires both U and V to be computed. Thin U and V are enough, there is no need for full U or V. 221*bf2c3715SXin Li * 222*bf2c3715SXin Li * \note SVD solving is implicitly least-squares. Thus, this method serves both purposes of exact solving and least-squares solving. 223*bf2c3715SXin Li * In other words, the returned solution is guaranteed to minimize the Euclidean norm \f$ \Vert A x - b \Vert \f$. 224*bf2c3715SXin Li */ 225*bf2c3715SXin Li template<typename Rhs> 226*bf2c3715SXin Li inline const Solve<Derived, Rhs> 227*bf2c3715SXin Li solve(const MatrixBase<Rhs>& b) const; 228*bf2c3715SXin Li #endif 229*bf2c3715SXin Li 230*bf2c3715SXin Li 231*bf2c3715SXin Li /** \brief Reports whether previous computation was successful. 232*bf2c3715SXin Li * 233*bf2c3715SXin Li * \returns \c Success if computation was successful. 234*bf2c3715SXin Li */ 235*bf2c3715SXin Li EIGEN_DEVICE_FUNC 236*bf2c3715SXin Li ComputationInfo info() const 237*bf2c3715SXin Li { 238*bf2c3715SXin Li eigen_assert(m_isInitialized && "SVD is not initialized."); 239*bf2c3715SXin Li return m_info; 240*bf2c3715SXin Li } 241*bf2c3715SXin Li 242*bf2c3715SXin Li #ifndef EIGEN_PARSED_BY_DOXYGEN 243*bf2c3715SXin Li template<typename RhsType, typename DstType> 244*bf2c3715SXin Li void _solve_impl(const RhsType &rhs, DstType &dst) const; 245*bf2c3715SXin Li 246*bf2c3715SXin Li template<bool Conjugate, typename RhsType, typename DstType> 247*bf2c3715SXin Li void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const; 248*bf2c3715SXin Li #endif 249*bf2c3715SXin Li 250*bf2c3715SXin Li protected: 251*bf2c3715SXin Li 252*bf2c3715SXin Li static void check_template_parameters() 253*bf2c3715SXin Li { 254*bf2c3715SXin Li EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar); 255*bf2c3715SXin Li } 256*bf2c3715SXin Li 257*bf2c3715SXin Li void _check_compute_assertions() const { 258*bf2c3715SXin Li eigen_assert(m_isInitialized && "SVD is not initialized."); 259*bf2c3715SXin Li } 260*bf2c3715SXin Li 261*bf2c3715SXin Li template<bool Transpose_, typename Rhs> 262*bf2c3715SXin Li void _check_solve_assertion(const Rhs& b) const { 263*bf2c3715SXin Li EIGEN_ONLY_USED_FOR_DEBUG(b); 264*bf2c3715SXin Li _check_compute_assertions(); 265*bf2c3715SXin Li eigen_assert(computeU() && computeV() && "SVDBase::solve(): Both unitaries U and V are required to be computed (thin unitaries suffice)."); 266*bf2c3715SXin Li eigen_assert((Transpose_?cols():rows())==b.rows() && "SVDBase::solve(): invalid number of rows of the right hand side matrix b"); 267*bf2c3715SXin Li } 268*bf2c3715SXin Li 269*bf2c3715SXin Li // return true if already allocated 270*bf2c3715SXin Li bool allocate(Index rows, Index cols, unsigned int computationOptions) ; 271*bf2c3715SXin Li 272*bf2c3715SXin Li MatrixUType m_matrixU; 273*bf2c3715SXin Li MatrixVType m_matrixV; 274*bf2c3715SXin Li SingularValuesType m_singularValues; 275*bf2c3715SXin Li ComputationInfo m_info; 276*bf2c3715SXin Li bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold; 277*bf2c3715SXin Li bool m_computeFullU, m_computeThinU; 278*bf2c3715SXin Li bool m_computeFullV, m_computeThinV; 279*bf2c3715SXin Li unsigned int m_computationOptions; 280*bf2c3715SXin Li Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize; 281*bf2c3715SXin Li RealScalar m_prescribedThreshold; 282*bf2c3715SXin Li 283*bf2c3715SXin Li /** \brief Default Constructor. 284*bf2c3715SXin Li * 285*bf2c3715SXin Li * Default constructor of SVDBase 286*bf2c3715SXin Li */ 287*bf2c3715SXin Li SVDBase() 288*bf2c3715SXin Li : m_info(Success), 289*bf2c3715SXin Li m_isInitialized(false), 290*bf2c3715SXin Li m_isAllocated(false), 291*bf2c3715SXin Li m_usePrescribedThreshold(false), 292*bf2c3715SXin Li m_computeFullU(false), 293*bf2c3715SXin Li m_computeThinU(false), 294*bf2c3715SXin Li m_computeFullV(false), 295*bf2c3715SXin Li m_computeThinV(false), 296*bf2c3715SXin Li m_computationOptions(0), 297*bf2c3715SXin Li m_rows(-1), m_cols(-1), m_diagSize(0) 298*bf2c3715SXin Li { 299*bf2c3715SXin Li check_template_parameters(); 300*bf2c3715SXin Li } 301*bf2c3715SXin Li 302*bf2c3715SXin Li 303*bf2c3715SXin Li }; 304*bf2c3715SXin Li 305*bf2c3715SXin Li #ifndef EIGEN_PARSED_BY_DOXYGEN 306*bf2c3715SXin Li template<typename Derived> 307*bf2c3715SXin Li template<typename RhsType, typename DstType> 308*bf2c3715SXin Li void SVDBase<Derived>::_solve_impl(const RhsType &rhs, DstType &dst) const 309*bf2c3715SXin Li { 310*bf2c3715SXin Li // A = U S V^* 311*bf2c3715SXin Li // So A^{-1} = V S^{-1} U^* 312*bf2c3715SXin Li 313*bf2c3715SXin Li Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp; 314*bf2c3715SXin Li Index l_rank = rank(); 315*bf2c3715SXin Li tmp.noalias() = m_matrixU.leftCols(l_rank).adjoint() * rhs; 316*bf2c3715SXin Li tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp; 317*bf2c3715SXin Li dst = m_matrixV.leftCols(l_rank) * tmp; 318*bf2c3715SXin Li } 319*bf2c3715SXin Li 320*bf2c3715SXin Li template<typename Derived> 321*bf2c3715SXin Li template<bool Conjugate, typename RhsType, typename DstType> 322*bf2c3715SXin Li void SVDBase<Derived>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const 323*bf2c3715SXin Li { 324*bf2c3715SXin Li // A = U S V^* 325*bf2c3715SXin Li // So A^{-*} = U S^{-1} V^* 326*bf2c3715SXin Li // And A^{-T} = U_conj S^{-1} V^T 327*bf2c3715SXin Li Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp; 328*bf2c3715SXin Li Index l_rank = rank(); 329*bf2c3715SXin Li 330*bf2c3715SXin Li tmp.noalias() = m_matrixV.leftCols(l_rank).transpose().template conjugateIf<Conjugate>() * rhs; 331*bf2c3715SXin Li tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp; 332*bf2c3715SXin Li dst = m_matrixU.template conjugateIf<!Conjugate>().leftCols(l_rank) * tmp; 333*bf2c3715SXin Li } 334*bf2c3715SXin Li #endif 335*bf2c3715SXin Li 336*bf2c3715SXin Li template<typename MatrixType> 337*bf2c3715SXin Li bool SVDBase<MatrixType>::allocate(Index rows, Index cols, unsigned int computationOptions) 338*bf2c3715SXin Li { 339*bf2c3715SXin Li eigen_assert(rows >= 0 && cols >= 0); 340*bf2c3715SXin Li 341*bf2c3715SXin Li if (m_isAllocated && 342*bf2c3715SXin Li rows == m_rows && 343*bf2c3715SXin Li cols == m_cols && 344*bf2c3715SXin Li computationOptions == m_computationOptions) 345*bf2c3715SXin Li { 346*bf2c3715SXin Li return true; 347*bf2c3715SXin Li } 348*bf2c3715SXin Li 349*bf2c3715SXin Li m_rows = rows; 350*bf2c3715SXin Li m_cols = cols; 351*bf2c3715SXin Li m_info = Success; 352*bf2c3715SXin Li m_isInitialized = false; 353*bf2c3715SXin Li m_isAllocated = true; 354*bf2c3715SXin Li m_computationOptions = computationOptions; 355*bf2c3715SXin Li m_computeFullU = (computationOptions & ComputeFullU) != 0; 356*bf2c3715SXin Li m_computeThinU = (computationOptions & ComputeThinU) != 0; 357*bf2c3715SXin Li m_computeFullV = (computationOptions & ComputeFullV) != 0; 358*bf2c3715SXin Li m_computeThinV = (computationOptions & ComputeThinV) != 0; 359*bf2c3715SXin Li eigen_assert(!(m_computeFullU && m_computeThinU) && "SVDBase: you can't ask for both full and thin U"); 360*bf2c3715SXin Li eigen_assert(!(m_computeFullV && m_computeThinV) && "SVDBase: you can't ask for both full and thin V"); 361*bf2c3715SXin Li eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) && 362*bf2c3715SXin Li "SVDBase: thin U and V are only available when your matrix has a dynamic number of columns."); 363*bf2c3715SXin Li 364*bf2c3715SXin Li m_diagSize = (std::min)(m_rows, m_cols); 365*bf2c3715SXin Li m_singularValues.resize(m_diagSize); 366*bf2c3715SXin Li if(RowsAtCompileTime==Dynamic) 367*bf2c3715SXin Li m_matrixU.resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0); 368*bf2c3715SXin Li if(ColsAtCompileTime==Dynamic) 369*bf2c3715SXin Li m_matrixV.resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0); 370*bf2c3715SXin Li 371*bf2c3715SXin Li return false; 372*bf2c3715SXin Li } 373*bf2c3715SXin Li 374*bf2c3715SXin Li }// end namespace 375*bf2c3715SXin Li 376*bf2c3715SXin Li #endif // EIGEN_SVDBASE_H 377