JacobiSVD.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2009-2010 Benoit Jacob <jacob.benoit.1@gmail.com>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_JACOBISVD_H
00011 #define EIGEN_JACOBISVD_H
00012 
00013 namespace Eigen { 
00014 
00015 namespace internal {
00016 // forward declaration (needed by ICC)
00017 // the empty body is required by MSVC
00018 template<typename MatrixType, int QRPreconditioner,
00019          bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
00020 struct svd_precondition_2x2_block_to_be_real {};
00021 
00022 /*** QR preconditioners (R-SVD)
00023  ***
00024  *** Their role is to reduce the problem of computing the SVD to the case of a square matrix.
00025  *** This approach, known as R-SVD, is an optimization for rectangular-enough matrices, and is a requirement for
00026  *** JacobiSVD which by itself is only able to work on square matrices.
00027  ***/
00028 
00029 enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
00030 
00031 template<typename MatrixType, int QRPreconditioner, int Case>
00032 struct qr_preconditioner_should_do_anything
00033 {
00034   enum { a = MatrixType::RowsAtCompileTime != Dynamic &&
00035              MatrixType::ColsAtCompileTime != Dynamic &&
00036              MatrixType::ColsAtCompileTime <= MatrixType::RowsAtCompileTime,
00037          b = MatrixType::RowsAtCompileTime != Dynamic &&
00038              MatrixType::ColsAtCompileTime != Dynamic &&
00039              MatrixType::RowsAtCompileTime <= MatrixType::ColsAtCompileTime,
00040          ret = !( (QRPreconditioner == NoQRPreconditioner) ||
00041                   (Case == PreconditionIfMoreColsThanRows && bool(a)) ||
00042                   (Case == PreconditionIfMoreRowsThanCols && bool(b)) )
00043   };
00044 };
00045 
00046 template<typename MatrixType, int QRPreconditioner, int Case,
00047          bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret
00048 > struct qr_preconditioner_impl {};
00049 
00050 template<typename MatrixType, int QRPreconditioner, int Case>
00051 class qr_preconditioner_impl<MatrixType, QRPreconditioner, Case, false>
00052 {
00053 public:
00054   typedef typename MatrixType::Index Index;
00055   void allocate(const JacobiSVD<MatrixType, QRPreconditioner>&) {}
00056   bool run(JacobiSVD<MatrixType, QRPreconditioner>&, const MatrixType&)
00057   {
00058     return false;
00059   }
00060 };
00061 
00062 /*** preconditioner using FullPivHouseholderQR ***/
00063 
00064 template<typename MatrixType>
00065 class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
00066 {
00067 public:
00068   typedef typename MatrixType::Index Index;
00069   typedef typename MatrixType::Scalar Scalar;
00070   enum
00071   {
00072     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00073     MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime
00074   };
00075   typedef Matrix<Scalar, 1, RowsAtCompileTime, RowMajor, 1, MaxRowsAtCompileTime> WorkspaceType;
00076 
00077   void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
00078   {
00079     if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
00080     {
00081       m_qr = FullPivHouseholderQR<MatrixType>(svd.rows(), svd.cols());
00082     }
00083     if (svd.m_computeFullU) m_workspace.resize(svd.rows());
00084   }
00085 
00086   bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00087   {
00088     if(matrix.rows() > matrix.cols())
00089     {
00090       m_qr.compute(matrix);
00091       svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00092       if(svd.m_computeFullU) m_qr.matrixQ().evalTo(svd.m_matrixU, m_workspace);
00093       if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
00094       return true;
00095     }
00096     return false;
00097   }
00098 private:
00099   FullPivHouseholderQR<MatrixType> m_qr;
00100   WorkspaceType m_workspace;
00101 };
00102 
00103 template<typename MatrixType>
00104 class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00105 {
00106 public:
00107   typedef typename MatrixType::Index Index;
00108   typedef typename MatrixType::Scalar Scalar;
00109   enum
00110   {
00111     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00112     ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00113     MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00114     MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00115     Options = MatrixType::Options
00116   };
00117   typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
00118           TransposeTypeWithSameStorageOrder;
00119 
00120   void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
00121   {
00122     if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
00123     {
00124       m_qr = FullPivHouseholderQR<TransposeTypeWithSameStorageOrder>(svd.cols(), svd.rows());
00125     }
00126     m_adjoint.resize(svd.cols(), svd.rows());
00127     if (svd.m_computeFullV) m_workspace.resize(svd.cols());
00128   }
00129 
00130   bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00131   {
00132     if(matrix.cols() > matrix.rows())
00133     {
00134       m_adjoint = matrix.adjoint();
00135       m_qr.compute(m_adjoint);
00136       svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00137       if(svd.m_computeFullV) m_qr.matrixQ().evalTo(svd.m_matrixV, m_workspace);
00138       if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
00139       return true;
00140     }
00141     else return false;
00142   }
00143 private:
00144   FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> m_qr;
00145   TransposeTypeWithSameStorageOrder m_adjoint;
00146   typename internal::plain_row_type<MatrixType>::type m_workspace;
00147 };
00148 
00149 /*** preconditioner using ColPivHouseholderQR ***/
00150 
00151 template<typename MatrixType>
00152 class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
00153 {
00154 public:
00155   typedef typename MatrixType::Index Index;
00156 
00157   void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
00158   {
00159     if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
00160     {
00161       m_qr = ColPivHouseholderQR<MatrixType>(svd.rows(), svd.cols());
00162     }
00163     if (svd.m_computeFullU) m_workspace.resize(svd.rows());
00164     else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
00165   }
00166 
00167   bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00168   {
00169     if(matrix.rows() > matrix.cols())
00170     {
00171       m_qr.compute(matrix);
00172       svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00173       if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
00174       else if(svd.m_computeThinU)
00175       {
00176         svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
00177         m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
00178       }
00179       if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
00180       return true;
00181     }
00182     return false;
00183   }
00184 
00185 private:
00186   ColPivHouseholderQR<MatrixType> m_qr;
00187   typename internal::plain_col_type<MatrixType>::type m_workspace;
00188 };
00189 
00190 template<typename MatrixType>
00191 class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00192 {
00193 public:
00194   typedef typename MatrixType::Index Index;
00195   typedef typename MatrixType::Scalar Scalar;
00196   enum
00197   {
00198     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00199     ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00200     MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00201     MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00202     Options = MatrixType::Options
00203   };
00204 
00205   typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
00206           TransposeTypeWithSameStorageOrder;
00207 
00208   void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
00209   {
00210     if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
00211     {
00212       m_qr = ColPivHouseholderQR<TransposeTypeWithSameStorageOrder>(svd.cols(), svd.rows());
00213     }
00214     if (svd.m_computeFullV) m_workspace.resize(svd.cols());
00215     else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
00216     m_adjoint.resize(svd.cols(), svd.rows());
00217   }
00218 
00219   bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00220   {
00221     if(matrix.cols() > matrix.rows())
00222     {
00223       m_adjoint = matrix.adjoint();
00224       m_qr.compute(m_adjoint);
00225 
00226       svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00227       if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
00228       else if(svd.m_computeThinV)
00229       {
00230         svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
00231         m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
00232       }
00233       if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
00234       return true;
00235     }
00236     else return false;
00237   }
00238 
00239 private:
00240   ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> m_qr;
00241   TransposeTypeWithSameStorageOrder m_adjoint;
00242   typename internal::plain_row_type<MatrixType>::type m_workspace;
00243 };
00244 
00245 /*** preconditioner using HouseholderQR ***/
00246 
00247 template<typename MatrixType>
00248 class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
00249 {
00250 public:
00251   typedef typename MatrixType::Index Index;
00252 
00253   void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
00254   {
00255     if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
00256     {
00257       m_qr = HouseholderQR<MatrixType>(svd.rows(), svd.cols());
00258     }
00259     if (svd.m_computeFullU) m_workspace.resize(svd.rows());
00260     else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
00261   }
00262 
00263   bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00264   {
00265     if(matrix.rows() > matrix.cols())
00266     {
00267       m_qr.compute(matrix);
00268       svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00269       if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
00270       else if(svd.m_computeThinU)
00271       {
00272         svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
00273         m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
00274       }
00275       if(svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
00276       return true;
00277     }
00278     return false;
00279   }
00280 private:
00281   HouseholderQR<MatrixType> m_qr;
00282   typename internal::plain_col_type<MatrixType>::type m_workspace;
00283 };
00284 
00285 template<typename MatrixType>
00286 class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00287 {
00288 public:
00289   typedef typename MatrixType::Index Index;
00290   typedef typename MatrixType::Scalar Scalar;
00291   enum
00292   {
00293     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00294     ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00295     MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00296     MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00297     Options = MatrixType::Options
00298   };
00299 
00300   typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
00301           TransposeTypeWithSameStorageOrder;
00302 
00303   void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
00304   {
00305     if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
00306     {
00307       m_qr = HouseholderQR<TransposeTypeWithSameStorageOrder>(svd.cols(), svd.rows());
00308     }
00309     if (svd.m_computeFullV) m_workspace.resize(svd.cols());
00310     else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
00311     m_adjoint.resize(svd.cols(), svd.rows());
00312   }
00313 
00314   bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00315   {
00316     if(matrix.cols() > matrix.rows())
00317     {
00318       m_adjoint = matrix.adjoint();
00319       m_qr.compute(m_adjoint);
00320 
00321       svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00322       if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
00323       else if(svd.m_computeThinV)
00324       {
00325         svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
00326         m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
00327       }
00328       if(svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
00329       return true;
00330     }
00331     else return false;
00332   }
00333 
00334 private:
00335   HouseholderQR<TransposeTypeWithSameStorageOrder> m_qr;
00336   TransposeTypeWithSameStorageOrder m_adjoint;
00337   typename internal::plain_row_type<MatrixType>::type m_workspace;
00338 };
00339 
00340 /*** 2x2 SVD implementation
00341  ***
00342  *** JacobiSVD consists in performing a series of 2x2 SVD subproblems
00343  ***/
00344 
00345 template<typename MatrixType, int QRPreconditioner>
00346 struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, false>
00347 {
00348   typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
00349   typedef typename SVD::Index Index;
00350   static void run(typename SVD::WorkMatrixType&, SVD&, Index, Index) {}
00351 };
00352 
00353 template<typename MatrixType, int QRPreconditioner>
00354 struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
00355 {
00356   typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
00357   typedef typename MatrixType::Scalar Scalar;
00358   typedef typename MatrixType::RealScalar RealScalar;
00359   typedef typename SVD::Index Index;
00360   static void run(typename SVD::WorkMatrixType& work_matrix, SVD& svd, Index p, Index q)
00361   {
00362     Scalar z;
00363     JacobiRotation<Scalar> rot;
00364     RealScalar n = sqrt(abs2(work_matrix.coeff(p,p)) + abs2(work_matrix.coeff(q,p)));
00365     if(n==0)
00366     {
00367       z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
00368       work_matrix.row(p) *= z;
00369       if(svd.computeU()) svd.m_matrixU.col(p) *= conj(z);
00370       z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
00371       work_matrix.row(q) *= z;
00372       if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
00373     }
00374     else
00375     {
00376       rot.c() = conj(work_matrix.coeff(p,p)) / n;
00377       rot.s() = work_matrix.coeff(q,p) / n;
00378       work_matrix.applyOnTheLeft(p,q,rot);
00379       if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
00380       if(work_matrix.coeff(p,q) != Scalar(0))
00381       {
00382         Scalar z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
00383         work_matrix.col(q) *= z;
00384         if(svd.computeV()) svd.m_matrixV.col(q) *= z;
00385       }
00386       if(work_matrix.coeff(q,q) != Scalar(0))
00387       {
00388         z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
00389         work_matrix.row(q) *= z;
00390         if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
00391       }
00392     }
00393   }
00394 };
00395 
00396 template<typename MatrixType, typename RealScalar, typename Index>
00397 void real_2x2_jacobi_svd(const MatrixType& matrix, Index p, Index q,
00398                             JacobiRotation<RealScalar> *j_left,
00399                             JacobiRotation<RealScalar> *j_right)
00400 {
00401   Matrix<RealScalar,2,2> m;
00402   m << real(matrix.coeff(p,p)), real(matrix.coeff(p,q)),
00403        real(matrix.coeff(q,p)), real(matrix.coeff(q,q));
00404   JacobiRotation<RealScalar> rot1;
00405   RealScalar t = m.coeff(0,0) + m.coeff(1,1);
00406   RealScalar d = m.coeff(1,0) - m.coeff(0,1);
00407   if(t == RealScalar(0))
00408   {
00409     rot1.c() = RealScalar(0);
00410     rot1.s() = d > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
00411   }
00412   else
00413   {
00414     RealScalar u = d / t;
00415     rot1.c() = RealScalar(1) / sqrt(RealScalar(1) + abs2(u));
00416     rot1.s() = rot1.c() * u;
00417   }
00418   m.applyOnTheLeft(0,1,rot1);
00419   j_right->makeJacobi(m,0,1);
00420   *j_left  = rot1 * j_right->transpose();
00421 }
00422 
00423 } // end namespace internal
00424 
00478 template<typename _MatrixType, int QRPreconditioner> class JacobiSVD
00479 {
00480   public:
00481 
00482     typedef _MatrixType MatrixType;
00483     typedef typename MatrixType::Scalar Scalar;
00484     typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
00485     typedef typename MatrixType::Index Index;
00486     enum {
00487       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00488       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00489       DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
00490       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00491       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00492       MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
00493       MatrixOptions = MatrixType::Options
00494     };
00495 
00496     typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
00497                    MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
00498             MatrixUType;
00499     typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime,
00500                    MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime>
00501             MatrixVType;
00502     typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
00503     typedef typename internal::plain_row_type<MatrixType>::type RowType;
00504     typedef typename internal::plain_col_type<MatrixType>::type ColType;
00505     typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
00506                    MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
00507             WorkMatrixType;
00508 
00514     JacobiSVD()
00515       : m_isInitialized(false),
00516         m_isAllocated(false),
00517         m_computationOptions(0),
00518         m_rows(-1), m_cols(-1)
00519     {}
00520 
00521 
00528     JacobiSVD(Index rows, Index cols, unsigned int computationOptions = 0)
00529       : m_isInitialized(false),
00530         m_isAllocated(false),
00531         m_computationOptions(0),
00532         m_rows(-1), m_cols(-1)
00533     {
00534       allocate(rows, cols, computationOptions);
00535     }
00536 
00547     JacobiSVD(const MatrixType& matrix, unsigned int computationOptions = 0)
00548       : m_isInitialized(false),
00549         m_isAllocated(false),
00550         m_computationOptions(0),
00551         m_rows(-1), m_cols(-1)
00552     {
00553       compute(matrix, computationOptions);
00554     }
00555 
00566     JacobiSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
00567 
00574     JacobiSVD& compute(const MatrixType& matrix)
00575     {
00576       return compute(matrix, m_computationOptions);
00577     }
00578 
00588     const MatrixUType& matrixU() const
00589     {
00590       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00591       eigen_assert(computeU() && "This JacobiSVD decomposition didn't compute U. Did you ask for it?");
00592       return m_matrixU;
00593     }
00594 
00604     const MatrixVType& matrixV() const
00605     {
00606       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00607       eigen_assert(computeV() && "This JacobiSVD decomposition didn't compute V. Did you ask for it?");
00608       return m_matrixV;
00609     }
00610 
00616     const SingularValuesType& singularValues() const
00617     {
00618       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00619       return m_singularValues;
00620     }
00621 
00623     inline bool computeU() const { return m_computeFullU || m_computeThinU; }
00625     inline bool computeV() const { return m_computeFullV || m_computeThinV; }
00626 
00636     template<typename Rhs>
00637     inline const internal::solve_retval<JacobiSVD, Rhs>
00638     solve(const MatrixBase<Rhs>& b) const
00639     {
00640       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00641       eigen_assert(computeU() && computeV() && "JacobiSVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
00642       return internal::solve_retval<JacobiSVD, Rhs>(*this, b.derived());
00643     }
00644 
00646     Index nonzeroSingularValues() const
00647     {
00648       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00649       return m_nonzeroSingularValues;
00650     }
00651 
00652     inline Index rows() const { return m_rows; }
00653     inline Index cols() const { return m_cols; }
00654 
00655   private:
00656     void allocate(Index rows, Index cols, unsigned int computationOptions);
00657 
00658   protected:
00659     MatrixUType m_matrixU;
00660     MatrixVType m_matrixV;
00661     SingularValuesType m_singularValues;
00662     WorkMatrixType m_workMatrix;
00663     bool m_isInitialized, m_isAllocated;
00664     bool m_computeFullU, m_computeThinU;
00665     bool m_computeFullV, m_computeThinV;
00666     unsigned int m_computationOptions;
00667     Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
00668 
00669     template<typename __MatrixType, int _QRPreconditioner, bool _IsComplex>
00670     friend struct internal::svd_precondition_2x2_block_to_be_real;
00671     template<typename __MatrixType, int _QRPreconditioner, int _Case, bool _DoAnything>
00672     friend struct internal::qr_preconditioner_impl;
00673 
00674     internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols;
00675     internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows;
00676 };
00677 
00678 template<typename MatrixType, int QRPreconditioner>
00679 void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols, unsigned int computationOptions)
00680 {
00681   eigen_assert(rows >= 0 && cols >= 0);
00682 
00683   if (m_isAllocated &&
00684       rows == m_rows &&
00685       cols == m_cols &&
00686       computationOptions == m_computationOptions)
00687   {
00688     return;
00689   }
00690 
00691   m_rows = rows;
00692   m_cols = cols;
00693   m_isInitialized = false;
00694   m_isAllocated = true;
00695   m_computationOptions = computationOptions;
00696   m_computeFullU = (computationOptions & ComputeFullU) != 0;
00697   m_computeThinU = (computationOptions & ComputeThinU) != 0;
00698   m_computeFullV = (computationOptions & ComputeFullV) != 0;
00699   m_computeThinV = (computationOptions & ComputeThinV) != 0;
00700   eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U");
00701   eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V");
00702   eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
00703               "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
00704   if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
00705   {
00706       eigen_assert(!(m_computeThinU || m_computeThinV) &&
00707               "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
00708               "Use the ColPivHouseholderQR preconditioner instead.");
00709   }
00710   m_diagSize = (std::min)(m_rows, m_cols);
00711   m_singularValues.resize(m_diagSize);
00712   m_matrixU.resize(m_rows, m_computeFullU ? m_rows
00713                           : m_computeThinU ? m_diagSize
00714                           : 0);
00715   m_matrixV.resize(m_cols, m_computeFullV ? m_cols
00716                           : m_computeThinV ? m_diagSize
00717                           : 0);
00718   m_workMatrix.resize(m_diagSize, m_diagSize);
00719   
00720   if(m_cols>m_rows) m_qr_precond_morecols.allocate(*this);
00721   if(m_rows>m_cols) m_qr_precond_morerows.allocate(*this);
00722 }
00723 
00724 template<typename MatrixType, int QRPreconditioner>
00725 JacobiSVD<MatrixType, QRPreconditioner>&
00726 JacobiSVD<MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)
00727 {
00728   allocate(matrix.rows(), matrix.cols(), computationOptions);
00729 
00730   // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
00731   // only worsening the precision of U and V as we accumulate more rotations
00732   const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
00733 
00734   // limit for very small denormal numbers to be considered zero in order to avoid infinite loops (see bug 286)
00735   const RealScalar considerAsZero = RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min();
00736 
00737   /*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */
00738 
00739   if(!m_qr_precond_morecols.run(*this, matrix) && !m_qr_precond_morerows.run(*this, matrix))
00740   {
00741     m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize);
00742     if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
00743     if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
00744     if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
00745     if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
00746   }
00747 
00748   /*** step 2. The main Jacobi SVD iteration. ***/
00749 
00750   bool finished = false;
00751   while(!finished)
00752   {
00753     finished = true;
00754 
00755     // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
00756 
00757     for(Index p = 1; p < m_diagSize; ++p)
00758     {
00759       for(Index q = 0; q < p; ++q)
00760       {
00761         // if this 2x2 sub-matrix is not diagonal already...
00762         // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
00763         // keep us iterating forever. Similarly, small denormal numbers are considered zero.
00764         using std::max;
00765         RealScalar threshold = (max)(considerAsZero, precision * (max)(internal::abs(m_workMatrix.coeff(p,p)),
00766                                                                        internal::abs(m_workMatrix.coeff(q,q))));
00767         if((max)(internal::abs(m_workMatrix.coeff(p,q)),internal::abs(m_workMatrix.coeff(q,p))) > threshold)
00768         {
00769           finished = false;
00770 
00771           // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
00772           internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q);
00773           JacobiRotation<RealScalar> j_left, j_right;
00774           internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
00775 
00776           // accumulate resulting Jacobi rotations
00777           m_workMatrix.applyOnTheLeft(p,q,j_left);
00778           if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
00779 
00780           m_workMatrix.applyOnTheRight(p,q,j_right);
00781           if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
00782         }
00783       }
00784     }
00785   }
00786 
00787   /*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/
00788 
00789   for(Index i = 0; i < m_diagSize; ++i)
00790   {
00791     RealScalar a = internal::abs(m_workMatrix.coeff(i,i));
00792     m_singularValues.coeffRef(i) = a;
00793     if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
00794   }
00795 
00796   /*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/
00797 
00798   m_nonzeroSingularValues = m_diagSize;
00799   for(Index i = 0; i < m_diagSize; i++)
00800   {
00801     Index pos;
00802     RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
00803     if(maxRemainingSingularValue == RealScalar(0))
00804     {
00805       m_nonzeroSingularValues = i;
00806       break;
00807     }
00808     if(pos)
00809     {
00810       pos += i;
00811       std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
00812       if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
00813       if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
00814     }
00815   }
00816 
00817   m_isInitialized = true;
00818   return *this;
00819 }
00820 
00821 namespace internal {
00822 template<typename _MatrixType, int QRPreconditioner, typename Rhs>
00823 struct solve_retval<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
00824   : solve_retval_base<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
00825 {
00826   typedef JacobiSVD<_MatrixType, QRPreconditioner> JacobiSVDType;
00827   EIGEN_MAKE_SOLVE_HELPERS(JacobiSVDType,Rhs)
00828 
00829   template<typename Dest> void evalTo(Dest& dst) const
00830   {
00831     eigen_assert(rhs().rows() == dec().rows());
00832 
00833     // A = U S V^*
00834     // So A^{-1} = V S^{-1} U^*
00835 
00836     Index diagSize = (std::min)(dec().rows(), dec().cols());
00837     typename JacobiSVDType::SingularValuesType invertedSingVals(diagSize);
00838 
00839     Index nonzeroSingVals = dec().nonzeroSingularValues();
00840     invertedSingVals.head(nonzeroSingVals) = dec().singularValues().head(nonzeroSingVals).array().inverse();
00841     invertedSingVals.tail(diagSize - nonzeroSingVals).setZero();
00842 
00843     dst = dec().matrixV().leftCols(diagSize)
00844         * invertedSingVals.asDiagonal()
00845         * dec().matrixU().leftCols(diagSize).adjoint()
00846         * rhs();
00847   }
00848 };
00849 } // end namespace internal
00850 
00858 template<typename Derived>
00859 JacobiSVD<typename MatrixBase<Derived>::PlainObject>
00860 MatrixBase<Derived>::jacobiSvd(unsigned int computationOptions) const
00861 {
00862   return JacobiSVD<PlainObject>(*this, computationOptions);
00863 }
00864 
00865 } // end namespace Eigen
00866 
00867 #endif // EIGEN_JACOBISVD_H


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:11:02