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 // Eigen is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 3 of the License, or (at your option) any later version.
00010 //
00011 // Alternatively, you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License as
00013 // published by the Free Software Foundation; either version 2 of
00014 // the License, or (at your option) any later version.
00015 //
00016 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00017 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00018 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00019 // GNU General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License and a copy of the GNU General Public License along with
00023 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00024 
00025 #ifndef EIGEN_JACOBISVD_H
00026 #define EIGEN_JACOBISVD_H
00027 
00028 namespace internal {
00029 // forward declaration (needed by ICC)
00030 // the empty body is required by MSVC
00031 template<typename MatrixType, int QRPreconditioner,
00032          bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
00033 struct svd_precondition_2x2_block_to_be_real {};
00034 
00035 /*** QR preconditioners (R-SVD)
00036  ***
00037  *** Their role is to reduce the problem of computing the SVD to the case of a square matrix.
00038  *** This approach, known as R-SVD, is an optimization for rectangular-enough matrices, and is a requirement for
00039  *** JacobiSVD which by itself is only able to work on square matrices.
00040  ***/
00041 
00042 enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
00043 
00044 template<typename MatrixType, int QRPreconditioner, int Case>
00045 struct qr_preconditioner_should_do_anything
00046 {
00047   enum { a = MatrixType::RowsAtCompileTime != Dynamic &&
00048              MatrixType::ColsAtCompileTime != Dynamic &&
00049              MatrixType::ColsAtCompileTime <= MatrixType::RowsAtCompileTime,
00050          b = MatrixType::RowsAtCompileTime != Dynamic &&
00051              MatrixType::ColsAtCompileTime != Dynamic &&
00052              MatrixType::RowsAtCompileTime <= MatrixType::ColsAtCompileTime,
00053          ret = !( (QRPreconditioner == NoQRPreconditioner) ||
00054                   (Case == PreconditionIfMoreColsThanRows && bool(a)) ||
00055                   (Case == PreconditionIfMoreRowsThanCols && bool(b)) )
00056   };
00057 };
00058 
00059 template<typename MatrixType, int QRPreconditioner, int Case,
00060          bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret
00061 > struct qr_preconditioner_impl {};
00062 
00063 template<typename MatrixType, int QRPreconditioner, int Case>
00064 struct qr_preconditioner_impl<MatrixType, QRPreconditioner, Case, false>
00065 {
00066   static bool run(JacobiSVD<MatrixType, QRPreconditioner>&, const MatrixType&)
00067   {
00068     return false;
00069   }
00070 };
00071 
00072 /*** preconditioner using FullPivHouseholderQR ***/
00073 
00074 template<typename MatrixType>
00075 struct qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
00076 {
00077   static bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00078   {
00079     if(matrix.rows() > matrix.cols())
00080     {
00081       FullPivHouseholderQR<MatrixType> qr(matrix);
00082       svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00083       if(svd.m_computeFullU) svd.m_matrixU = qr.matrixQ();
00084       if(svd.computeV()) svd.m_matrixV = qr.colsPermutation();
00085       return true;
00086     }
00087     return false;
00088   }
00089 };
00090 
00091 template<typename MatrixType>
00092 struct qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00093 {
00094   static bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00095   {
00096     if(matrix.cols() > matrix.rows())
00097     {
00098       typedef Matrix<typename MatrixType::Scalar, MatrixType::ColsAtCompileTime, MatrixType::RowsAtCompileTime,
00099                      MatrixType::Options, MatrixType::MaxColsAtCompileTime, MatrixType::MaxRowsAtCompileTime>
00100               TransposeTypeWithSameStorageOrder;
00101       FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> qr(matrix.adjoint());
00102       svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00103       if(svd.m_computeFullV) svd.m_matrixV = qr.matrixQ();
00104       if(svd.computeU()) svd.m_matrixU = qr.colsPermutation();
00105       return true;
00106     }
00107     else return false;
00108   }
00109 };
00110 
00111 /*** preconditioner using ColPivHouseholderQR ***/
00112 
00113 template<typename MatrixType>
00114 struct qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
00115 {
00116   static bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00117   {
00118     if(matrix.rows() > matrix.cols())
00119     {
00120       ColPivHouseholderQR<MatrixType> qr(matrix);
00121       svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00122       if(svd.m_computeFullU) svd.m_matrixU = qr.householderQ();
00123       else if(svd.m_computeThinU) {
00124         svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
00125         qr.householderQ().applyThisOnTheLeft(svd.m_matrixU);
00126       }
00127       if(svd.computeV()) svd.m_matrixV = qr.colsPermutation();
00128       return true;
00129     }
00130     return false;
00131   }
00132 };
00133 
00134 template<typename MatrixType>
00135 struct qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00136 {
00137   static bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00138   {
00139     if(matrix.cols() > matrix.rows())
00140     {
00141       typedef Matrix<typename MatrixType::Scalar, MatrixType::ColsAtCompileTime, MatrixType::RowsAtCompileTime,
00142                      MatrixType::Options, MatrixType::MaxColsAtCompileTime, MatrixType::MaxRowsAtCompileTime>
00143               TransposeTypeWithSameStorageOrder;
00144       ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> qr(matrix.adjoint());
00145       svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00146       if(svd.m_computeFullV) svd.m_matrixV = qr.householderQ();
00147       else if(svd.m_computeThinV) {
00148         svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
00149         qr.householderQ().applyThisOnTheLeft(svd.m_matrixV);
00150       }
00151       if(svd.computeU()) svd.m_matrixU = qr.colsPermutation();
00152       return true;
00153     }
00154     else return false;
00155   }
00156 };
00157 
00158 /*** preconditioner using HouseholderQR ***/
00159 
00160 template<typename MatrixType>
00161 struct qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
00162 {
00163   static bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00164   {
00165     if(matrix.rows() > matrix.cols())
00166     {
00167       HouseholderQR<MatrixType> qr(matrix);
00168       svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00169       if(svd.m_computeFullU) svd.m_matrixU = qr.householderQ();
00170       else if(svd.m_computeThinU) {
00171         svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
00172         qr.householderQ().applyThisOnTheLeft(svd.m_matrixU);
00173       }
00174       if(svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
00175       return true;
00176     }
00177     return false;
00178   }
00179 };
00180 
00181 template<typename MatrixType>
00182 struct qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00183 {
00184   static bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00185   {
00186     if(matrix.cols() > matrix.rows())
00187     {
00188       typedef Matrix<typename MatrixType::Scalar, MatrixType::ColsAtCompileTime, MatrixType::RowsAtCompileTime,
00189                      MatrixType::Options, MatrixType::MaxColsAtCompileTime, MatrixType::MaxRowsAtCompileTime>
00190               TransposeTypeWithSameStorageOrder;
00191       HouseholderQR<TransposeTypeWithSameStorageOrder> qr(matrix.adjoint());
00192       svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00193       if(svd.m_computeFullV) svd.m_matrixV = qr.householderQ();
00194       else if(svd.m_computeThinV) {
00195         svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
00196         qr.householderQ().applyThisOnTheLeft(svd.m_matrixV);
00197       }
00198       if(svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
00199       return true;
00200     }
00201     else return false;
00202   }
00203 };
00204 
00205 /*** 2x2 SVD implementation
00206  ***
00207  *** JacobiSVD consists in performing a series of 2x2 SVD subproblems
00208  ***/
00209 
00210 template<typename MatrixType, int QRPreconditioner>
00211 struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, false>
00212 {
00213   typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
00214   typedef typename SVD::Index Index;
00215   static void run(typename SVD::WorkMatrixType&, SVD&, Index, Index) {}
00216 };
00217 
00218 template<typename MatrixType, int QRPreconditioner>
00219 struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
00220 {
00221   typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
00222   typedef typename MatrixType::Scalar Scalar;
00223   typedef typename MatrixType::RealScalar RealScalar;
00224   typedef typename SVD::Index Index;
00225   static void run(typename SVD::WorkMatrixType& work_matrix, SVD& svd, Index p, Index q)
00226   {
00227     Scalar z;
00228     JacobiRotation<Scalar> rot;
00229     RealScalar n = sqrt(abs2(work_matrix.coeff(p,p)) + abs2(work_matrix.coeff(q,p)));
00230     if(n==0)
00231     {
00232       z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
00233       work_matrix.row(p) *= z;
00234       if(svd.computeU()) svd.m_matrixU.col(p) *= conj(z);
00235       z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
00236       work_matrix.row(q) *= z;
00237       if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
00238     }
00239     else
00240     {
00241       rot.c() = conj(work_matrix.coeff(p,p)) / n;
00242       rot.s() = work_matrix.coeff(q,p) / n;
00243       work_matrix.applyOnTheLeft(p,q,rot);
00244       if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
00245       if(work_matrix.coeff(p,q) != Scalar(0))
00246       {
00247         Scalar z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
00248         work_matrix.col(q) *= z;
00249         if(svd.computeV()) svd.m_matrixV.col(q) *= z;
00250       }
00251       if(work_matrix.coeff(q,q) != Scalar(0))
00252       {
00253         z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
00254         work_matrix.row(q) *= z;
00255         if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
00256       }
00257     }
00258   }
00259 };
00260 
00261 template<typename MatrixType, typename RealScalar, typename Index>
00262 void real_2x2_jacobi_svd(const MatrixType& matrix, Index p, Index q,
00263                             JacobiRotation<RealScalar> *j_left,
00264                             JacobiRotation<RealScalar> *j_right)
00265 {
00266   Matrix<RealScalar,2,2> m;
00267   m << real(matrix.coeff(p,p)), real(matrix.coeff(p,q)),
00268        real(matrix.coeff(q,p)), real(matrix.coeff(q,q));
00269   JacobiRotation<RealScalar> rot1;
00270   RealScalar t = m.coeff(0,0) + m.coeff(1,1);
00271   RealScalar d = m.coeff(1,0) - m.coeff(0,1);
00272   if(t == RealScalar(0))
00273   {
00274     rot1.c() = RealScalar(0);
00275     rot1.s() = d > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
00276   }
00277   else
00278   {
00279     RealScalar u = d / t;
00280     rot1.c() = RealScalar(1) / sqrt(RealScalar(1) + abs2(u));
00281     rot1.s() = rot1.c() * u;
00282   }
00283   m.applyOnTheLeft(0,1,rot1);
00284   j_right->makeJacobi(m,0,1);
00285   *j_left  = rot1 * j_right->transpose();
00286 }
00287 
00288 } // end namespace internal
00289 
00343 template<typename _MatrixType, int QRPreconditioner> class JacobiSVD
00344 {
00345   public:
00346 
00347     typedef _MatrixType MatrixType;
00348     typedef typename MatrixType::Scalar Scalar;
00349     typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
00350     typedef typename MatrixType::Index Index;
00351     enum {
00352       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00353       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00354       DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
00355       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00356       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00357       MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
00358       MatrixOptions = MatrixType::Options
00359     };
00360 
00361     typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
00362                    MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
00363             MatrixUType;
00364     typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime,
00365                    MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime>
00366             MatrixVType;
00367     typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
00368     typedef typename internal::plain_row_type<MatrixType>::type RowType;
00369     typedef typename internal::plain_col_type<MatrixType>::type ColType;
00370     typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
00371                    MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
00372             WorkMatrixType;
00373 
00379     JacobiSVD()
00380       : m_isInitialized(false),
00381         m_isAllocated(false),
00382         m_computationOptions(0),
00383         m_rows(-1), m_cols(-1)
00384     {}
00385 
00386 
00393     JacobiSVD(Index rows, Index cols, unsigned int computationOptions = 0)
00394       : m_isInitialized(false),
00395         m_isAllocated(false),
00396         m_computationOptions(0),
00397         m_rows(-1), m_cols(-1)
00398     {
00399       allocate(rows, cols, computationOptions);
00400     }
00401 
00412     JacobiSVD(const MatrixType& matrix, unsigned int computationOptions = 0)
00413       : m_isInitialized(false),
00414         m_isAllocated(false),
00415         m_computationOptions(0),
00416         m_rows(-1), m_cols(-1)
00417     {
00418       compute(matrix, computationOptions);
00419     }
00420 
00431     JacobiSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
00432 
00439     JacobiSVD& compute(const MatrixType& matrix)
00440     {
00441       return compute(matrix, m_computationOptions);
00442     }
00443 
00453     const MatrixUType& matrixU() const
00454     {
00455       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00456       eigen_assert(computeU() && "This JacobiSVD decomposition didn't compute U. Did you ask for it?");
00457       return m_matrixU;
00458     }
00459 
00469     const MatrixVType& matrixV() const
00470     {
00471       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00472       eigen_assert(computeV() && "This JacobiSVD decomposition didn't compute V. Did you ask for it?");
00473       return m_matrixV;
00474     }
00475 
00481     const SingularValuesType& singularValues() const
00482     {
00483       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00484       return m_singularValues;
00485     }
00486 
00488     inline bool computeU() const { return m_computeFullU || m_computeThinU; }
00490     inline bool computeV() const { return m_computeFullV || m_computeThinV; }
00491 
00501     template<typename Rhs>
00502     inline const internal::solve_retval<JacobiSVD, Rhs>
00503     solve(const MatrixBase<Rhs>& b) const
00504     {
00505       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00506       eigen_assert(computeU() && computeV() && "JacobiSVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
00507       return internal::solve_retval<JacobiSVD, Rhs>(*this, b.derived());
00508     }
00509 
00511     Index nonzeroSingularValues() const
00512     {
00513       eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00514       return m_nonzeroSingularValues;
00515     }
00516 
00517     inline Index rows() const { return m_rows; }
00518     inline Index cols() const { return m_cols; }
00519 
00520   private:
00521     void allocate(Index rows, Index cols, unsigned int computationOptions);
00522 
00523   protected:
00524     MatrixUType m_matrixU;
00525     MatrixVType m_matrixV;
00526     SingularValuesType m_singularValues;
00527     WorkMatrixType m_workMatrix;
00528     bool m_isInitialized, m_isAllocated;
00529     bool m_computeFullU, m_computeThinU;
00530     bool m_computeFullV, m_computeThinV;
00531     unsigned int m_computationOptions;
00532     Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
00533 
00534     template<typename __MatrixType, int _QRPreconditioner, bool _IsComplex>
00535     friend struct internal::svd_precondition_2x2_block_to_be_real;
00536     template<typename __MatrixType, int _QRPreconditioner, int _Case, bool _DoAnything>
00537     friend struct internal::qr_preconditioner_impl;
00538 };
00539 
00540 template<typename MatrixType, int QRPreconditioner>
00541 void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols, unsigned int computationOptions)
00542 {
00543   eigen_assert(rows >= 0 && cols >= 0);
00544 
00545   if (m_isAllocated &&
00546       rows == m_rows &&
00547       cols == m_cols &&
00548       computationOptions == m_computationOptions)
00549   {
00550     return;
00551   }
00552 
00553   m_rows = rows;
00554   m_cols = cols;
00555   m_isInitialized = false;
00556   m_isAllocated = true;
00557   m_computationOptions = computationOptions;
00558   m_computeFullU = (computationOptions & ComputeFullU) != 0;
00559   m_computeThinU = (computationOptions & ComputeThinU) != 0;
00560   m_computeFullV = (computationOptions & ComputeFullV) != 0;
00561   m_computeThinV = (computationOptions & ComputeThinV) != 0;
00562   eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U");
00563   eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V");
00564   eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
00565               "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
00566   if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
00567   {
00568       eigen_assert(!(m_computeThinU || m_computeThinV) &&
00569               "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
00570               "Use the ColPivHouseholderQR preconditioner instead.");
00571   }
00572   m_diagSize = (std::min)(m_rows, m_cols);
00573   m_singularValues.resize(m_diagSize);
00574   m_matrixU.resize(m_rows, m_computeFullU ? m_rows
00575                           : m_computeThinU ? m_diagSize
00576                           : 0);
00577   m_matrixV.resize(m_cols, m_computeFullV ? m_cols
00578                           : m_computeThinV ? m_diagSize
00579                           : 0);
00580   m_workMatrix.resize(m_diagSize, m_diagSize);
00581 }
00582 
00583 template<typename MatrixType, int QRPreconditioner>
00584 JacobiSVD<MatrixType, QRPreconditioner>&
00585 JacobiSVD<MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)
00586 {
00587   allocate(matrix.rows(), matrix.cols(), computationOptions);
00588 
00589   // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
00590   // only worsening the precision of U and V as we accumulate more rotations
00591   const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
00592 
00593   /*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */
00594 
00595   if(!internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows>::run(*this, matrix)
00596   && !internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols>::run(*this, matrix))
00597   {
00598     m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize);
00599     if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
00600     if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
00601     if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
00602     if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
00603   }
00604 
00605   /*** step 2. The main Jacobi SVD iteration. ***/
00606 
00607   bool finished = false;
00608   while(!finished)
00609   {
00610     finished = true;
00611 
00612     // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
00613 
00614     for(Index p = 1; p < m_diagSize; ++p)
00615     {
00616       for(Index q = 0; q < p; ++q)
00617       {
00618         // if this 2x2 sub-matrix is not diagonal already...
00619         // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
00620         // keep us iterating forever.
00621         using std::max;
00622         if((max)(internal::abs(m_workMatrix.coeff(p,q)),internal::abs(m_workMatrix.coeff(q,p)))
00623             > (max)(internal::abs(m_workMatrix.coeff(p,p)),internal::abs(m_workMatrix.coeff(q,q)))*precision)
00624         {
00625           finished = false;
00626 
00627           // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
00628           internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q);
00629           JacobiRotation<RealScalar> j_left, j_right;
00630           internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
00631 
00632           // accumulate resulting Jacobi rotations
00633           m_workMatrix.applyOnTheLeft(p,q,j_left);
00634           if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
00635 
00636           m_workMatrix.applyOnTheRight(p,q,j_right);
00637           if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
00638         }
00639       }
00640     }
00641   }
00642 
00643   /*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/
00644 
00645   for(Index i = 0; i < m_diagSize; ++i)
00646   {
00647     RealScalar a = internal::abs(m_workMatrix.coeff(i,i));
00648     m_singularValues.coeffRef(i) = a;
00649     if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
00650   }
00651 
00652   /*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/
00653 
00654   m_nonzeroSingularValues = m_diagSize;
00655   for(Index i = 0; i < m_diagSize; i++)
00656   {
00657     Index pos;
00658     RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
00659     if(maxRemainingSingularValue == RealScalar(0))
00660     {
00661       m_nonzeroSingularValues = i;
00662       break;
00663     }
00664     if(pos)
00665     {
00666       pos += i;
00667       std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
00668       if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
00669       if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
00670     }
00671   }
00672 
00673   m_isInitialized = true;
00674   return *this;
00675 }
00676 
00677 namespace internal {
00678 template<typename _MatrixType, int QRPreconditioner, typename Rhs>
00679 struct solve_retval<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
00680   : solve_retval_base<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
00681 {
00682   typedef JacobiSVD<_MatrixType, QRPreconditioner> JacobiSVDType;
00683   EIGEN_MAKE_SOLVE_HELPERS(JacobiSVDType,Rhs)
00684 
00685   template<typename Dest> void evalTo(Dest& dst) const
00686   {
00687     eigen_assert(rhs().rows() == dec().rows());
00688 
00689     // A = U S V^*
00690     // So A^{-1} = V S^{-1} U^*
00691 
00692     Index diagSize = (std::min)(dec().rows(), dec().cols());
00693     typename JacobiSVDType::SingularValuesType invertedSingVals(diagSize);
00694 
00695     Index nonzeroSingVals = dec().nonzeroSingularValues();
00696     invertedSingVals.head(nonzeroSingVals) = dec().singularValues().head(nonzeroSingVals).array().inverse();
00697     invertedSingVals.tail(diagSize - nonzeroSingVals).setZero();
00698 
00699     dst = dec().matrixV().leftCols(diagSize)
00700         * invertedSingVals.asDiagonal()
00701         * dec().matrixU().leftCols(diagSize).adjoint()
00702         * rhs();
00703   }
00704 };
00705 } // end namespace internal
00706 
00707 template<typename Derived>
00708 JacobiSVD<typename MatrixBase<Derived>::PlainObject>
00709 MatrixBase<Derived>::jacobiSvd(unsigned int computationOptions) const
00710 {
00711   return JacobiSVD<PlainObject>(*this, computationOptions);
00712 }
00713 
00714 
00715 
00716 #endif // EIGEN_JACOBISVD_H


libicr
Author(s): Robert Krug
autogenerated on Mon Jan 6 2014 11:32:53