Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/Eigen/src/SVD/SVDBase.h was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2009-2010 Benoit Jacob <jacob.benoit.1@gmail.com>
0005 // Copyright (C) 2014 Gael Guennebaud <gael.guennebaud@inria.fr>
0006 //
0007 // Copyright (C) 2013 Gauthier Brun <brun.gauthier@gmail.com>
0008 // Copyright (C) 2013 Nicolas Carre <nicolas.carre@ensimag.fr>
0009 // Copyright (C) 2013 Jean Ceccato <jean.ceccato@ensimag.fr>
0010 // Copyright (C) 2013 Pierre Zoppitelli <pierre.zoppitelli@ensimag.fr>
0011 //
0012 // This Source Code Form is subject to the terms of the Mozilla
0013 // Public License v. 2.0. If a copy of the MPL was not distributed
0014 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0015 
0016 #ifndef EIGEN_SVDBASE_H
0017 #define EIGEN_SVDBASE_H
0018 
0019 namespace Eigen {
0020 
0021 namespace internal {
0022 template<typename Derived> struct traits<SVDBase<Derived> >
0023  : traits<Derived>
0024 {
0025   typedef MatrixXpr XprKind;
0026   typedef SolverStorage StorageKind;
0027   typedef int StorageIndex;
0028   enum { Flags = 0 };
0029 };
0030 }
0031 
0032 /** \ingroup SVD_Module
0033  *
0034  *
0035  * \class SVDBase
0036  *
0037  * \brief Base class of SVD algorithms
0038  *
0039  * \tparam Derived the type of the actual SVD decomposition
0040  *
0041  * SVD decomposition consists in decomposing any n-by-p matrix \a A as a product
0042  *   \f[ A = U S V^* \f]
0043  * 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;
0044  * 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
0045  * and right \em singular \em vectors of \a A respectively.
0046  *
0047  * Singular values are always sorted in decreasing order.
0048  *
0049  * 
0050  * 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
0051  * 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
0052  * 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,
0053  * 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.
0054  * 
0055  * The status of the computation can be retrived using the \a info() method. Unless \a info() returns \a Success, the results should be not
0056  * considered well defined.
0057  *  
0058  * 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
0059  * terminate in finite (and reasonable) time.
0060  * \sa class BDCSVD, class JacobiSVD
0061  */
0062 template<typename Derived> class SVDBase
0063  : public SolverBase<SVDBase<Derived> >
0064 {
0065 public: 
0066    
0067   template<typename Derived_>
0068   friend struct internal::solve_assertion;
0069 
0070   typedef typename internal::traits<Derived>::MatrixType MatrixType;
0071   typedef typename MatrixType::Scalar Scalar;
0072   typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
0073   typedef typename Eigen::internal::traits<SVDBase>::StorageIndex StorageIndex;
0074   typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
0075   enum {
0076     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
0077     ColsAtCompileTime = MatrixType::ColsAtCompileTime,
0078     DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
0079     MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
0080     MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
0081     MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
0082     MatrixOptions = MatrixType::Options
0083   };
0084 
0085   typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixUType;
0086   typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime> MatrixVType;
0087   typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
0088   
0089   Derived& derived() { return *static_cast<Derived*>(this); }
0090   const Derived& derived() const { return *static_cast<const Derived*>(this); }
0091 
0092   /** \returns the \a U matrix.
0093    *
0094    * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p,
0095    * 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.
0096    *
0097    * The \a m first columns of \a U are the left singular vectors of the matrix being decomposed.
0098    *
0099    * This method asserts that you asked for \a U to be computed.
0100    */
0101   const MatrixUType& matrixU() const
0102   {
0103     _check_compute_assertions();
0104     eigen_assert(computeU() && "This SVD decomposition didn't compute U. Did you ask for it?");
0105     return m_matrixU;
0106   }
0107 
0108   /** \returns the \a V matrix.
0109    *
0110    * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p,
0111    * 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.
0112    *
0113    * The \a m first columns of \a V are the right singular vectors of the matrix being decomposed.
0114    *
0115    * This method asserts that you asked for \a V to be computed.
0116    */
0117   const MatrixVType& matrixV() const
0118   {
0119     _check_compute_assertions();
0120     eigen_assert(computeV() && "This SVD decomposition didn't compute V. Did you ask for it?");
0121     return m_matrixV;
0122   }
0123 
0124   /** \returns the vector of singular values.
0125    *
0126    * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p, the
0127    * returned vector has size \a m.  Singular values are always sorted in decreasing order.
0128    */
0129   const SingularValuesType& singularValues() const
0130   {
0131     _check_compute_assertions();
0132     return m_singularValues;
0133   }
0134 
0135   /** \returns the number of singular values that are not exactly 0 */
0136   Index nonzeroSingularValues() const
0137   {
0138     _check_compute_assertions();
0139     return m_nonzeroSingularValues;
0140   }
0141   
0142   /** \returns the rank of the matrix of which \c *this is the SVD.
0143     *
0144     * \note This method has to determine which singular values should be considered nonzero.
0145     *       For that, it uses the threshold value that you can control by calling
0146     *       setThreshold(const RealScalar&).
0147     */
0148   inline Index rank() const
0149   {
0150     using std::abs;
0151     _check_compute_assertions();
0152     if(m_singularValues.size()==0) return 0;
0153     RealScalar premultiplied_threshold = numext::maxi<RealScalar>(m_singularValues.coeff(0) * threshold(), (std::numeric_limits<RealScalar>::min)());
0154     Index i = m_nonzeroSingularValues-1;
0155     while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
0156     return i+1;
0157   }
0158   
0159   /** Allows to prescribe a threshold to be used by certain methods, such as rank() and solve(),
0160     * which need to determine when singular values are to be considered nonzero.
0161     * This is not used for the SVD decomposition itself.
0162     *
0163     * When it needs to get the threshold value, Eigen calls threshold().
0164     * The default is \c NumTraits<Scalar>::epsilon()
0165     *
0166     * \param threshold The new value to use as the threshold.
0167     *
0168     * A singular value will be considered nonzero if its value is strictly greater than
0169     *  \f$ \vert singular value \vert \leqslant threshold \times \vert max singular value \vert \f$.
0170     *
0171     * If you want to come back to the default behavior, call setThreshold(Default_t)
0172     */
0173   Derived& setThreshold(const RealScalar& threshold)
0174   {
0175     m_usePrescribedThreshold = true;
0176     m_prescribedThreshold = threshold;
0177     return derived();
0178   }
0179 
0180   /** Allows to come back to the default behavior, letting Eigen use its default formula for
0181     * determining the threshold.
0182     *
0183     * You should pass the special object Eigen::Default as parameter here.
0184     * \code svd.setThreshold(Eigen::Default); \endcode
0185     *
0186     * See the documentation of setThreshold(const RealScalar&).
0187     */
0188   Derived& setThreshold(Default_t)
0189   {
0190     m_usePrescribedThreshold = false;
0191     return derived();
0192   }
0193 
0194   /** Returns the threshold that will be used by certain methods such as rank().
0195     *
0196     * See the documentation of setThreshold(const RealScalar&).
0197     */
0198   RealScalar threshold() const
0199   {
0200     eigen_assert(m_isInitialized || m_usePrescribedThreshold);
0201     // this temporary is needed to workaround a MSVC issue
0202     Index diagSize = (std::max<Index>)(1,m_diagSize);
0203     return m_usePrescribedThreshold ? m_prescribedThreshold
0204                                     : RealScalar(diagSize)*NumTraits<Scalar>::epsilon();
0205   }
0206 
0207   /** \returns true if \a U (full or thin) is asked for in this SVD decomposition */
0208   inline bool computeU() const { return m_computeFullU || m_computeThinU; }
0209   /** \returns true if \a V (full or thin) is asked for in this SVD decomposition */
0210   inline bool computeV() const { return m_computeFullV || m_computeThinV; }
0211 
0212   inline Index rows() const { return m_rows; }
0213   inline Index cols() const { return m_cols; }
0214   
0215   #ifdef EIGEN_PARSED_BY_DOXYGEN
0216   /** \returns a (least squares) solution of \f$ A x = b \f$ using the current SVD decomposition of A.
0217     *
0218     * \param b the right-hand-side of the equation to solve.
0219     *
0220     * \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.
0221     *
0222     * \note SVD solving is implicitly least-squares. Thus, this method serves both purposes of exact solving and least-squares solving.
0223     * In other words, the returned solution is guaranteed to minimize the Euclidean norm \f$ \Vert A x - b \Vert \f$.
0224     */
0225   template<typename Rhs>
0226   inline const Solve<Derived, Rhs>
0227   solve(const MatrixBase<Rhs>& b) const;
0228   #endif
0229 
0230 
0231   /** \brief Reports whether previous computation was successful.
0232    *
0233    * \returns \c Success if computation was successful.
0234    */
0235   EIGEN_DEVICE_FUNC
0236   ComputationInfo info() const
0237   {
0238     eigen_assert(m_isInitialized && "SVD is not initialized.");
0239     return m_info;
0240   }
0241 
0242   #ifndef EIGEN_PARSED_BY_DOXYGEN
0243   template<typename RhsType, typename DstType>
0244   void _solve_impl(const RhsType &rhs, DstType &dst) const;
0245 
0246   template<bool Conjugate, typename RhsType, typename DstType>
0247   void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const;
0248   #endif
0249 
0250 protected:
0251 
0252   static void check_template_parameters()
0253   {
0254     EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
0255   }
0256 
0257   void _check_compute_assertions() const {
0258     eigen_assert(m_isInitialized && "SVD is not initialized.");
0259   }
0260 
0261   template<bool Transpose_, typename Rhs>
0262   void _check_solve_assertion(const Rhs& b) const {
0263       EIGEN_ONLY_USED_FOR_DEBUG(b);
0264       _check_compute_assertions();
0265       eigen_assert(computeU() && computeV() && "SVDBase::solve(): Both unitaries U and V are required to be computed (thin unitaries suffice).");
0266       eigen_assert((Transpose_?cols():rows())==b.rows() && "SVDBase::solve(): invalid number of rows of the right hand side matrix b");
0267   }
0268 
0269   // return true if already allocated
0270   bool allocate(Index rows, Index cols, unsigned int computationOptions) ;
0271 
0272   MatrixUType m_matrixU;
0273   MatrixVType m_matrixV;
0274   SingularValuesType m_singularValues;
0275   ComputationInfo m_info;
0276   bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
0277   bool m_computeFullU, m_computeThinU;
0278   bool m_computeFullV, m_computeThinV;
0279   unsigned int m_computationOptions;
0280   Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
0281   RealScalar m_prescribedThreshold;
0282 
0283   /** \brief Default Constructor.
0284    *
0285    * Default constructor of SVDBase
0286    */
0287   SVDBase()
0288     : m_info(Success),
0289       m_isInitialized(false),
0290       m_isAllocated(false),
0291       m_usePrescribedThreshold(false),
0292       m_computeFullU(false),
0293       m_computeThinU(false),
0294       m_computeFullV(false),
0295       m_computeThinV(false),
0296       m_computationOptions(0),
0297       m_rows(-1), m_cols(-1), m_diagSize(0)
0298   {
0299     check_template_parameters();
0300   }
0301 
0302 
0303 };
0304 
0305 #ifndef EIGEN_PARSED_BY_DOXYGEN
0306 template<typename Derived>
0307 template<typename RhsType, typename DstType>
0308 void SVDBase<Derived>::_solve_impl(const RhsType &rhs, DstType &dst) const
0309 {
0310   // A = U S V^*
0311   // So A^{-1} = V S^{-1} U^*
0312 
0313   Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp;
0314   Index l_rank = rank();
0315   tmp.noalias() =  m_matrixU.leftCols(l_rank).adjoint() * rhs;
0316   tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
0317   dst = m_matrixV.leftCols(l_rank) * tmp;
0318 }
0319 
0320 template<typename Derived>
0321 template<bool Conjugate, typename RhsType, typename DstType>
0322 void SVDBase<Derived>::_solve_impl_transposed(const RhsType &rhs, DstType &dst) const
0323 {
0324   // A = U S V^*
0325   // So  A^{-*} = U S^{-1} V^*
0326   // And A^{-T} = U_conj S^{-1} V^T
0327   Matrix<typename RhsType::Scalar, Dynamic, RhsType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, RhsType::MaxColsAtCompileTime> tmp;
0328   Index l_rank = rank();
0329 
0330   tmp.noalias() =  m_matrixV.leftCols(l_rank).transpose().template conjugateIf<Conjugate>() * rhs;
0331   tmp = m_singularValues.head(l_rank).asDiagonal().inverse() * tmp;
0332   dst = m_matrixU.template conjugateIf<!Conjugate>().leftCols(l_rank) * tmp;
0333 }
0334 #endif
0335 
0336 template<typename MatrixType>
0337 bool SVDBase<MatrixType>::allocate(Index rows, Index cols, unsigned int computationOptions)
0338 {
0339   eigen_assert(rows >= 0 && cols >= 0);
0340 
0341   if (m_isAllocated &&
0342       rows == m_rows &&
0343       cols == m_cols &&
0344       computationOptions == m_computationOptions)
0345   {
0346     return true;
0347   }
0348 
0349   m_rows = rows;
0350   m_cols = cols;
0351   m_info = Success;
0352   m_isInitialized = false;
0353   m_isAllocated = true;
0354   m_computationOptions = computationOptions;
0355   m_computeFullU = (computationOptions & ComputeFullU) != 0;
0356   m_computeThinU = (computationOptions & ComputeThinU) != 0;
0357   m_computeFullV = (computationOptions & ComputeFullV) != 0;
0358   m_computeThinV = (computationOptions & ComputeThinV) != 0;
0359   eigen_assert(!(m_computeFullU && m_computeThinU) && "SVDBase: you can't ask for both full and thin U");
0360   eigen_assert(!(m_computeFullV && m_computeThinV) && "SVDBase: you can't ask for both full and thin V");
0361   eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
0362            "SVDBase: thin U and V are only available when your matrix has a dynamic number of columns.");
0363 
0364   m_diagSize = (std::min)(m_rows, m_cols);
0365   m_singularValues.resize(m_diagSize);
0366   if(RowsAtCompileTime==Dynamic)
0367     m_matrixU.resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0);
0368   if(ColsAtCompileTime==Dynamic)
0369     m_matrixV.resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0);
0370 
0371   return false;
0372 }
0373 
0374 }// end namespace
0375 
0376 #endif // EIGEN_SVDBASE_H