00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #ifndef EIGEN_JACOBISVD_H
00011 #define EIGEN_JACOBISVD_H
00012
00013 namespace Eigen {
00014
00015 namespace internal {
00016
00017
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
00023
00024
00025
00026
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
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.~QRType();
00082 ::new (&m_qr) QRType(svd.rows(), svd.cols());
00083 }
00084 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
00085 }
00086
00087 bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00088 {
00089 if(matrix.rows() > matrix.cols())
00090 {
00091 m_qr.compute(matrix);
00092 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00093 if(svd.m_computeFullU) m_qr.matrixQ().evalTo(svd.m_matrixU, m_workspace);
00094 if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
00095 return true;
00096 }
00097 return false;
00098 }
00099 private:
00100 typedef FullPivHouseholderQR<MatrixType> QRType;
00101 QRType m_qr;
00102 WorkspaceType m_workspace;
00103 };
00104
00105 template<typename MatrixType>
00106 class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00107 {
00108 public:
00109 typedef typename MatrixType::Index Index;
00110 typedef typename MatrixType::Scalar Scalar;
00111 enum
00112 {
00113 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00114 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00115 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00116 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00117 Options = MatrixType::Options
00118 };
00119 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
00120 TransposeTypeWithSameStorageOrder;
00121
00122 void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
00123 {
00124 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
00125 {
00126 m_qr.~QRType();
00127 ::new (&m_qr) QRType(svd.cols(), svd.rows());
00128 }
00129 m_adjoint.resize(svd.cols(), svd.rows());
00130 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
00131 }
00132
00133 bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00134 {
00135 if(matrix.cols() > matrix.rows())
00136 {
00137 m_adjoint = matrix.adjoint();
00138 m_qr.compute(m_adjoint);
00139 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00140 if(svd.m_computeFullV) m_qr.matrixQ().evalTo(svd.m_matrixV, m_workspace);
00141 if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
00142 return true;
00143 }
00144 else return false;
00145 }
00146 private:
00147 typedef FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
00148 QRType m_qr;
00149 TransposeTypeWithSameStorageOrder m_adjoint;
00150 typename internal::plain_row_type<MatrixType>::type m_workspace;
00151 };
00152
00153
00154
00155 template<typename MatrixType>
00156 class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
00157 {
00158 public:
00159 typedef typename MatrixType::Index Index;
00160
00161 void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
00162 {
00163 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
00164 {
00165 m_qr.~QRType();
00166 ::new (&m_qr) QRType(svd.rows(), svd.cols());
00167 }
00168 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
00169 else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
00170 }
00171
00172 bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00173 {
00174 if(matrix.rows() > matrix.cols())
00175 {
00176 m_qr.compute(matrix);
00177 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00178 if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
00179 else if(svd.m_computeThinU)
00180 {
00181 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
00182 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
00183 }
00184 if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
00185 return true;
00186 }
00187 return false;
00188 }
00189
00190 private:
00191 typedef ColPivHouseholderQR<MatrixType> QRType;
00192 QRType m_qr;
00193 typename internal::plain_col_type<MatrixType>::type m_workspace;
00194 };
00195
00196 template<typename MatrixType>
00197 class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00198 {
00199 public:
00200 typedef typename MatrixType::Index Index;
00201 typedef typename MatrixType::Scalar Scalar;
00202 enum
00203 {
00204 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00205 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00206 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00207 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00208 Options = MatrixType::Options
00209 };
00210
00211 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
00212 TransposeTypeWithSameStorageOrder;
00213
00214 void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
00215 {
00216 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
00217 {
00218 m_qr.~QRType();
00219 ::new (&m_qr) QRType(svd.cols(), svd.rows());
00220 }
00221 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
00222 else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
00223 m_adjoint.resize(svd.cols(), svd.rows());
00224 }
00225
00226 bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00227 {
00228 if(matrix.cols() > matrix.rows())
00229 {
00230 m_adjoint = matrix.adjoint();
00231 m_qr.compute(m_adjoint);
00232
00233 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00234 if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
00235 else if(svd.m_computeThinV)
00236 {
00237 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
00238 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
00239 }
00240 if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
00241 return true;
00242 }
00243 else return false;
00244 }
00245
00246 private:
00247 typedef ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
00248 QRType m_qr;
00249 TransposeTypeWithSameStorageOrder m_adjoint;
00250 typename internal::plain_row_type<MatrixType>::type m_workspace;
00251 };
00252
00253
00254
00255 template<typename MatrixType>
00256 class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
00257 {
00258 public:
00259 typedef typename MatrixType::Index Index;
00260
00261 void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
00262 {
00263 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
00264 {
00265 m_qr.~QRType();
00266 ::new (&m_qr) QRType(svd.rows(), svd.cols());
00267 }
00268 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
00269 else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
00270 }
00271
00272 bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00273 {
00274 if(matrix.rows() > matrix.cols())
00275 {
00276 m_qr.compute(matrix);
00277 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
00278 if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
00279 else if(svd.m_computeThinU)
00280 {
00281 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
00282 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
00283 }
00284 if(svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
00285 return true;
00286 }
00287 return false;
00288 }
00289 private:
00290 typedef HouseholderQR<MatrixType> QRType;
00291 QRType m_qr;
00292 typename internal::plain_col_type<MatrixType>::type m_workspace;
00293 };
00294
00295 template<typename MatrixType>
00296 class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
00297 {
00298 public:
00299 typedef typename MatrixType::Index Index;
00300 typedef typename MatrixType::Scalar Scalar;
00301 enum
00302 {
00303 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00304 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00305 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00306 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00307 Options = MatrixType::Options
00308 };
00309
00310 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
00311 TransposeTypeWithSameStorageOrder;
00312
00313 void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
00314 {
00315 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
00316 {
00317 m_qr.~QRType();
00318 ::new (&m_qr) QRType(svd.cols(), svd.rows());
00319 }
00320 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
00321 else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
00322 m_adjoint.resize(svd.cols(), svd.rows());
00323 }
00324
00325 bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
00326 {
00327 if(matrix.cols() > matrix.rows())
00328 {
00329 m_adjoint = matrix.adjoint();
00330 m_qr.compute(m_adjoint);
00331
00332 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
00333 if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
00334 else if(svd.m_computeThinV)
00335 {
00336 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
00337 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
00338 }
00339 if(svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
00340 return true;
00341 }
00342 else return false;
00343 }
00344
00345 private:
00346 typedef HouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
00347 QRType m_qr;
00348 TransposeTypeWithSameStorageOrder m_adjoint;
00349 typename internal::plain_row_type<MatrixType>::type m_workspace;
00350 };
00351
00352
00353
00354
00355
00356
00357 template<typename MatrixType, int QRPreconditioner>
00358 struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, false>
00359 {
00360 typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
00361 typedef typename SVD::Index Index;
00362 static void run(typename SVD::WorkMatrixType&, SVD&, Index, Index) {}
00363 };
00364
00365 template<typename MatrixType, int QRPreconditioner>
00366 struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
00367 {
00368 typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
00369 typedef typename MatrixType::Scalar Scalar;
00370 typedef typename MatrixType::RealScalar RealScalar;
00371 typedef typename SVD::Index Index;
00372 static void run(typename SVD::WorkMatrixType& work_matrix, SVD& svd, Index p, Index q)
00373 {
00374 using std::sqrt;
00375 Scalar z;
00376 JacobiRotation<Scalar> rot;
00377 RealScalar n = sqrt(numext::abs2(work_matrix.coeff(p,p)) + numext::abs2(work_matrix.coeff(q,p)));
00378 if(n==0)
00379 {
00380 z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
00381 work_matrix.row(p) *= z;
00382 if(svd.computeU()) svd.m_matrixU.col(p) *= conj(z);
00383 z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
00384 work_matrix.row(q) *= z;
00385 if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
00386 }
00387 else
00388 {
00389 rot.c() = conj(work_matrix.coeff(p,p)) / n;
00390 rot.s() = work_matrix.coeff(q,p) / n;
00391 work_matrix.applyOnTheLeft(p,q,rot);
00392 if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
00393 if(work_matrix.coeff(p,q) != Scalar(0))
00394 {
00395 Scalar z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
00396 work_matrix.col(q) *= z;
00397 if(svd.computeV()) svd.m_matrixV.col(q) *= z;
00398 }
00399 if(work_matrix.coeff(q,q) != Scalar(0))
00400 {
00401 z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
00402 work_matrix.row(q) *= z;
00403 if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
00404 }
00405 }
00406 }
00407 };
00408
00409 template<typename MatrixType, typename RealScalar, typename Index>
00410 void real_2x2_jacobi_svd(const MatrixType& matrix, Index p, Index q,
00411 JacobiRotation<RealScalar> *j_left,
00412 JacobiRotation<RealScalar> *j_right)
00413 {
00414 using std::sqrt;
00415 Matrix<RealScalar,2,2> m;
00416 m << numext::real(matrix.coeff(p,p)), numext::real(matrix.coeff(p,q)),
00417 numext::real(matrix.coeff(q,p)), numext::real(matrix.coeff(q,q));
00418 JacobiRotation<RealScalar> rot1;
00419 RealScalar t = m.coeff(0,0) + m.coeff(1,1);
00420 RealScalar d = m.coeff(1,0) - m.coeff(0,1);
00421 if(t == RealScalar(0))
00422 {
00423 rot1.c() = RealScalar(0);
00424 rot1.s() = d > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
00425 }
00426 else
00427 {
00428 RealScalar u = d / t;
00429 rot1.c() = RealScalar(1) / sqrt(RealScalar(1) + numext::abs2(u));
00430 rot1.s() = rot1.c() * u;
00431 }
00432 m.applyOnTheLeft(0,1,rot1);
00433 j_right->makeJacobi(m,0,1);
00434 *j_left = rot1 * j_right->transpose();
00435 }
00436
00437 }
00438
00492 template<typename _MatrixType, int QRPreconditioner> class JacobiSVD
00493 {
00494 public:
00495
00496 typedef _MatrixType MatrixType;
00497 typedef typename MatrixType::Scalar Scalar;
00498 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
00499 typedef typename MatrixType::Index Index;
00500 enum {
00501 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00502 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00503 DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
00504 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00505 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00506 MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
00507 MatrixOptions = MatrixType::Options
00508 };
00509
00510 typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
00511 MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
00512 MatrixUType;
00513 typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime,
00514 MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime>
00515 MatrixVType;
00516 typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
00517 typedef typename internal::plain_row_type<MatrixType>::type RowType;
00518 typedef typename internal::plain_col_type<MatrixType>::type ColType;
00519 typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
00520 MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
00521 WorkMatrixType;
00522
00528 JacobiSVD()
00529 : m_isInitialized(false),
00530 m_isAllocated(false),
00531 m_computationOptions(0),
00532 m_rows(-1), m_cols(-1)
00533 {}
00534
00535
00542 JacobiSVD(Index rows, Index cols, unsigned int computationOptions = 0)
00543 : m_isInitialized(false),
00544 m_isAllocated(false),
00545 m_computationOptions(0),
00546 m_rows(-1), m_cols(-1)
00547 {
00548 allocate(rows, cols, computationOptions);
00549 }
00550
00561 JacobiSVD(const MatrixType& matrix, unsigned int computationOptions = 0)
00562 : m_isInitialized(false),
00563 m_isAllocated(false),
00564 m_computationOptions(0),
00565 m_rows(-1), m_cols(-1)
00566 {
00567 compute(matrix, computationOptions);
00568 }
00569
00580 JacobiSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
00581
00588 JacobiSVD& compute(const MatrixType& matrix)
00589 {
00590 return compute(matrix, m_computationOptions);
00591 }
00592
00602 const MatrixUType& matrixU() const
00603 {
00604 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00605 eigen_assert(computeU() && "This JacobiSVD decomposition didn't compute U. Did you ask for it?");
00606 return m_matrixU;
00607 }
00608
00618 const MatrixVType& matrixV() const
00619 {
00620 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00621 eigen_assert(computeV() && "This JacobiSVD decomposition didn't compute V. Did you ask for it?");
00622 return m_matrixV;
00623 }
00624
00630 const SingularValuesType& singularValues() const
00631 {
00632 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00633 return m_singularValues;
00634 }
00635
00637 inline bool computeU() const { return m_computeFullU || m_computeThinU; }
00639 inline bool computeV() const { return m_computeFullV || m_computeThinV; }
00640
00650 template<typename Rhs>
00651 inline const internal::solve_retval<JacobiSVD, Rhs>
00652 solve(const MatrixBase<Rhs>& b) const
00653 {
00654 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00655 eigen_assert(computeU() && computeV() && "JacobiSVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
00656 return internal::solve_retval<JacobiSVD, Rhs>(*this, b.derived());
00657 }
00658
00660 Index nonzeroSingularValues() const
00661 {
00662 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
00663 return m_nonzeroSingularValues;
00664 }
00665
00666 inline Index rows() const { return m_rows; }
00667 inline Index cols() const { return m_cols; }
00668
00669 private:
00670 void allocate(Index rows, Index cols, unsigned int computationOptions);
00671
00672 protected:
00673 MatrixUType m_matrixU;
00674 MatrixVType m_matrixV;
00675 SingularValuesType m_singularValues;
00676 WorkMatrixType m_workMatrix;
00677 bool m_isInitialized, m_isAllocated;
00678 bool m_computeFullU, m_computeThinU;
00679 bool m_computeFullV, m_computeThinV;
00680 unsigned int m_computationOptions;
00681 Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
00682
00683 template<typename __MatrixType, int _QRPreconditioner, bool _IsComplex>
00684 friend struct internal::svd_precondition_2x2_block_to_be_real;
00685 template<typename __MatrixType, int _QRPreconditioner, int _Case, bool _DoAnything>
00686 friend struct internal::qr_preconditioner_impl;
00687
00688 internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols;
00689 internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows;
00690 };
00691
00692 template<typename MatrixType, int QRPreconditioner>
00693 void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols, unsigned int computationOptions)
00694 {
00695 eigen_assert(rows >= 0 && cols >= 0);
00696
00697 if (m_isAllocated &&
00698 rows == m_rows &&
00699 cols == m_cols &&
00700 computationOptions == m_computationOptions)
00701 {
00702 return;
00703 }
00704
00705 m_rows = rows;
00706 m_cols = cols;
00707 m_isInitialized = false;
00708 m_isAllocated = true;
00709 m_computationOptions = computationOptions;
00710 m_computeFullU = (computationOptions & ComputeFullU) != 0;
00711 m_computeThinU = (computationOptions & ComputeThinU) != 0;
00712 m_computeFullV = (computationOptions & ComputeFullV) != 0;
00713 m_computeThinV = (computationOptions & ComputeThinV) != 0;
00714 eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U");
00715 eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V");
00716 eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
00717 "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
00718 if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
00719 {
00720 eigen_assert(!(m_computeThinU || m_computeThinV) &&
00721 "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
00722 "Use the ColPivHouseholderQR preconditioner instead.");
00723 }
00724 m_diagSize = (std::min)(m_rows, m_cols);
00725 m_singularValues.resize(m_diagSize);
00726 if(RowsAtCompileTime==Dynamic)
00727 m_matrixU.resize(m_rows, m_computeFullU ? m_rows
00728 : m_computeThinU ? m_diagSize
00729 : 0);
00730 if(ColsAtCompileTime==Dynamic)
00731 m_matrixV.resize(m_cols, m_computeFullV ? m_cols
00732 : m_computeThinV ? m_diagSize
00733 : 0);
00734 m_workMatrix.resize(m_diagSize, m_diagSize);
00735
00736 if(m_cols>m_rows) m_qr_precond_morecols.allocate(*this);
00737 if(m_rows>m_cols) m_qr_precond_morerows.allocate(*this);
00738 }
00739
00740 template<typename MatrixType, int QRPreconditioner>
00741 JacobiSVD<MatrixType, QRPreconditioner>&
00742 JacobiSVD<MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)
00743 {
00744 using std::abs;
00745 allocate(matrix.rows(), matrix.cols(), computationOptions);
00746
00747
00748
00749 const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
00750
00751
00752 const RealScalar considerAsZero = RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min();
00753
00754
00755
00756 if(!m_qr_precond_morecols.run(*this, matrix) && !m_qr_precond_morerows.run(*this, matrix))
00757 {
00758 m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize);
00759 if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
00760 if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
00761 if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
00762 if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
00763 }
00764
00765
00766
00767 bool finished = false;
00768 while(!finished)
00769 {
00770 finished = true;
00771
00772
00773
00774 for(Index p = 1; p < m_diagSize; ++p)
00775 {
00776 for(Index q = 0; q < p; ++q)
00777 {
00778
00779
00780
00781 using std::max;
00782 RealScalar threshold = (max)(considerAsZero, precision * (max)(abs(m_workMatrix.coeff(p,p)),
00783 abs(m_workMatrix.coeff(q,q))));
00784 if((max)(abs(m_workMatrix.coeff(p,q)),abs(m_workMatrix.coeff(q,p))) > threshold)
00785 {
00786 finished = false;
00787
00788
00789 internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q);
00790 JacobiRotation<RealScalar> j_left, j_right;
00791 internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
00792
00793
00794 m_workMatrix.applyOnTheLeft(p,q,j_left);
00795 if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
00796
00797 m_workMatrix.applyOnTheRight(p,q,j_right);
00798 if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
00799 }
00800 }
00801 }
00802 }
00803
00804
00805
00806 for(Index i = 0; i < m_diagSize; ++i)
00807 {
00808 RealScalar a = abs(m_workMatrix.coeff(i,i));
00809 m_singularValues.coeffRef(i) = a;
00810 if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
00811 }
00812
00813
00814
00815 m_nonzeroSingularValues = m_diagSize;
00816 for(Index i = 0; i < m_diagSize; i++)
00817 {
00818 Index pos;
00819 RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
00820 if(maxRemainingSingularValue == RealScalar(0))
00821 {
00822 m_nonzeroSingularValues = i;
00823 break;
00824 }
00825 if(pos)
00826 {
00827 pos += i;
00828 std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
00829 if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
00830 if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
00831 }
00832 }
00833
00834 m_isInitialized = true;
00835 return *this;
00836 }
00837
00838 namespace internal {
00839 template<typename _MatrixType, int QRPreconditioner, typename Rhs>
00840 struct solve_retval<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
00841 : solve_retval_base<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
00842 {
00843 typedef JacobiSVD<_MatrixType, QRPreconditioner> JacobiSVDType;
00844 EIGEN_MAKE_SOLVE_HELPERS(JacobiSVDType,Rhs)
00845
00846 template<typename Dest> void evalTo(Dest& dst) const
00847 {
00848 eigen_assert(rhs().rows() == dec().rows());
00849
00850
00851
00852
00853 Matrix<Scalar, Dynamic, Rhs::ColsAtCompileTime, 0, _MatrixType::MaxRowsAtCompileTime, Rhs::MaxColsAtCompileTime> tmp;
00854 Index nonzeroSingVals = dec().nonzeroSingularValues();
00855
00856 tmp.noalias() = dec().matrixU().leftCols(nonzeroSingVals).adjoint() * rhs();
00857 tmp = dec().singularValues().head(nonzeroSingVals).asDiagonal().inverse() * tmp;
00858 dst = dec().matrixV().leftCols(nonzeroSingVals) * tmp;
00859 }
00860 };
00861 }
00862
00870 template<typename Derived>
00871 JacobiSVD<typename MatrixBase<Derived>::PlainObject>
00872 MatrixBase<Derived>::jacobiSvd(unsigned int computationOptions) const
00873 {
00874 return JacobiSVD<PlainObject>(*this, computationOptions);
00875 }
00876
00877 }
00878
00879 #endif // EIGEN_JACOBISVD_H