1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2009-2010 Benoit Jacob <jacob.benoit.1@gmail.com>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_JACOBISVD_H
11#define EIGEN_JACOBISVD_H
12
13namespace Eigen {
14
15namespace internal {
16// forward declaration (needed by ICC)
17// the empty body is required by MSVC
18template<typename MatrixType, int QRPreconditioner,
19         bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
20struct svd_precondition_2x2_block_to_be_real {};
21
22/*** QR preconditioners (R-SVD)
23 ***
24 *** Their role is to reduce the problem of computing the SVD to the case of a square matrix.
25 *** This approach, known as R-SVD, is an optimization for rectangular-enough matrices, and is a requirement for
26 *** JacobiSVD which by itself is only able to work on square matrices.
27 ***/
28
29enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
30
31template<typename MatrixType, int QRPreconditioner, int Case>
32struct qr_preconditioner_should_do_anything
33{
34  enum { a = MatrixType::RowsAtCompileTime != Dynamic &&
35             MatrixType::ColsAtCompileTime != Dynamic &&
36             MatrixType::ColsAtCompileTime <= MatrixType::RowsAtCompileTime,
37         b = MatrixType::RowsAtCompileTime != Dynamic &&
38             MatrixType::ColsAtCompileTime != Dynamic &&
39             MatrixType::RowsAtCompileTime <= MatrixType::ColsAtCompileTime,
40         ret = !( (QRPreconditioner == NoQRPreconditioner) ||
41                  (Case == PreconditionIfMoreColsThanRows && bool(a)) ||
42                  (Case == PreconditionIfMoreRowsThanCols && bool(b)) )
43  };
44};
45
46template<typename MatrixType, int QRPreconditioner, int Case,
47         bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret
48> struct qr_preconditioner_impl {};
49
50template<typename MatrixType, int QRPreconditioner, int Case>
51class qr_preconditioner_impl<MatrixType, QRPreconditioner, Case, false>
52{
53public:
54  typedef typename MatrixType::Index Index;
55  void allocate(const JacobiSVD<MatrixType, QRPreconditioner>&) {}
56  bool run(JacobiSVD<MatrixType, QRPreconditioner>&, const MatrixType&)
57  {
58    return false;
59  }
60};
61
62/*** preconditioner using FullPivHouseholderQR ***/
63
64template<typename MatrixType>
65class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
66{
67public:
68  typedef typename MatrixType::Index Index;
69  typedef typename MatrixType::Scalar Scalar;
70  enum
71  {
72    RowsAtCompileTime = MatrixType::RowsAtCompileTime,
73    MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime
74  };
75  typedef Matrix<Scalar, 1, RowsAtCompileTime, RowMajor, 1, MaxRowsAtCompileTime> WorkspaceType;
76
77  void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
78  {
79    if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
80    {
81      m_qr.~QRType();
82      ::new (&m_qr) QRType(svd.rows(), svd.cols());
83    }
84    if (svd.m_computeFullU) m_workspace.resize(svd.rows());
85  }
86
87  bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
88  {
89    if(matrix.rows() > matrix.cols())
90    {
91      m_qr.compute(matrix);
92      svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
93      if(svd.m_computeFullU) m_qr.matrixQ().evalTo(svd.m_matrixU, m_workspace);
94      if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
95      return true;
96    }
97    return false;
98  }
99private:
100  typedef FullPivHouseholderQR<MatrixType> QRType;
101  QRType m_qr;
102  WorkspaceType m_workspace;
103};
104
105template<typename MatrixType>
106class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
107{
108public:
109  typedef typename MatrixType::Index Index;
110  typedef typename MatrixType::Scalar Scalar;
111  enum
112  {
113    RowsAtCompileTime = MatrixType::RowsAtCompileTime,
114    ColsAtCompileTime = MatrixType::ColsAtCompileTime,
115    MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
116    MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
117    Options = MatrixType::Options
118  };
119  typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
120          TransposeTypeWithSameStorageOrder;
121
122  void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
123  {
124    if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
125    {
126      m_qr.~QRType();
127      ::new (&m_qr) QRType(svd.cols(), svd.rows());
128    }
129    m_adjoint.resize(svd.cols(), svd.rows());
130    if (svd.m_computeFullV) m_workspace.resize(svd.cols());
131  }
132
133  bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
134  {
135    if(matrix.cols() > matrix.rows())
136    {
137      m_adjoint = matrix.adjoint();
138      m_qr.compute(m_adjoint);
139      svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
140      if(svd.m_computeFullV) m_qr.matrixQ().evalTo(svd.m_matrixV, m_workspace);
141      if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
142      return true;
143    }
144    else return false;
145  }
146private:
147  typedef FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
148  QRType m_qr;
149  TransposeTypeWithSameStorageOrder m_adjoint;
150  typename internal::plain_row_type<MatrixType>::type m_workspace;
151};
152
153/*** preconditioner using ColPivHouseholderQR ***/
154
155template<typename MatrixType>
156class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
157{
158public:
159  typedef typename MatrixType::Index Index;
160
161  void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
162  {
163    if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
164    {
165      m_qr.~QRType();
166      ::new (&m_qr) QRType(svd.rows(), svd.cols());
167    }
168    if (svd.m_computeFullU) m_workspace.resize(svd.rows());
169    else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
170  }
171
172  bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
173  {
174    if(matrix.rows() > matrix.cols())
175    {
176      m_qr.compute(matrix);
177      svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
178      if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
179      else if(svd.m_computeThinU)
180      {
181        svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
182        m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
183      }
184      if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
185      return true;
186    }
187    return false;
188  }
189
190private:
191  typedef ColPivHouseholderQR<MatrixType> QRType;
192  QRType m_qr;
193  typename internal::plain_col_type<MatrixType>::type m_workspace;
194};
195
196template<typename MatrixType>
197class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
198{
199public:
200  typedef typename MatrixType::Index Index;
201  typedef typename MatrixType::Scalar Scalar;
202  enum
203  {
204    RowsAtCompileTime = MatrixType::RowsAtCompileTime,
205    ColsAtCompileTime = MatrixType::ColsAtCompileTime,
206    MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
207    MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
208    Options = MatrixType::Options
209  };
210
211  typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
212          TransposeTypeWithSameStorageOrder;
213
214  void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
215  {
216    if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
217    {
218      m_qr.~QRType();
219      ::new (&m_qr) QRType(svd.cols(), svd.rows());
220    }
221    if (svd.m_computeFullV) m_workspace.resize(svd.cols());
222    else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
223    m_adjoint.resize(svd.cols(), svd.rows());
224  }
225
226  bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
227  {
228    if(matrix.cols() > matrix.rows())
229    {
230      m_adjoint = matrix.adjoint();
231      m_qr.compute(m_adjoint);
232
233      svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
234      if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
235      else if(svd.m_computeThinV)
236      {
237        svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
238        m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
239      }
240      if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
241      return true;
242    }
243    else return false;
244  }
245
246private:
247  typedef ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
248  QRType m_qr;
249  TransposeTypeWithSameStorageOrder m_adjoint;
250  typename internal::plain_row_type<MatrixType>::type m_workspace;
251};
252
253/*** preconditioner using HouseholderQR ***/
254
255template<typename MatrixType>
256class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
257{
258public:
259  typedef typename MatrixType::Index Index;
260
261  void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
262  {
263    if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
264    {
265      m_qr.~QRType();
266      ::new (&m_qr) QRType(svd.rows(), svd.cols());
267    }
268    if (svd.m_computeFullU) m_workspace.resize(svd.rows());
269    else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
270  }
271
272  bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
273  {
274    if(matrix.rows() > matrix.cols())
275    {
276      m_qr.compute(matrix);
277      svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
278      if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
279      else if(svd.m_computeThinU)
280      {
281        svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
282        m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
283      }
284      if(svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
285      return true;
286    }
287    return false;
288  }
289private:
290  typedef HouseholderQR<MatrixType> QRType;
291  QRType m_qr;
292  typename internal::plain_col_type<MatrixType>::type m_workspace;
293};
294
295template<typename MatrixType>
296class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
297{
298public:
299  typedef typename MatrixType::Index Index;
300  typedef typename MatrixType::Scalar Scalar;
301  enum
302  {
303    RowsAtCompileTime = MatrixType::RowsAtCompileTime,
304    ColsAtCompileTime = MatrixType::ColsAtCompileTime,
305    MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
306    MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
307    Options = MatrixType::Options
308  };
309
310  typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
311          TransposeTypeWithSameStorageOrder;
312
313  void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
314  {
315    if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
316    {
317      m_qr.~QRType();
318      ::new (&m_qr) QRType(svd.cols(), svd.rows());
319    }
320    if (svd.m_computeFullV) m_workspace.resize(svd.cols());
321    else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
322    m_adjoint.resize(svd.cols(), svd.rows());
323  }
324
325  bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
326  {
327    if(matrix.cols() > matrix.rows())
328    {
329      m_adjoint = matrix.adjoint();
330      m_qr.compute(m_adjoint);
331
332      svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
333      if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
334      else if(svd.m_computeThinV)
335      {
336        svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
337        m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
338      }
339      if(svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
340      return true;
341    }
342    else return false;
343  }
344
345private:
346  typedef HouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
347  QRType m_qr;
348  TransposeTypeWithSameStorageOrder m_adjoint;
349  typename internal::plain_row_type<MatrixType>::type m_workspace;
350};
351
352/*** 2x2 SVD implementation
353 ***
354 *** JacobiSVD consists in performing a series of 2x2 SVD subproblems
355 ***/
356
357template<typename MatrixType, int QRPreconditioner>
358struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, false>
359{
360  typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
361  typedef typename SVD::Index Index;
362  static void run(typename SVD::WorkMatrixType&, SVD&, Index, Index) {}
363};
364
365template<typename MatrixType, int QRPreconditioner>
366struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
367{
368  typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
369  typedef typename MatrixType::Scalar Scalar;
370  typedef typename MatrixType::RealScalar RealScalar;
371  typedef typename SVD::Index Index;
372  static void run(typename SVD::WorkMatrixType& work_matrix, SVD& svd, Index p, Index q)
373  {
374    using std::sqrt;
375    Scalar z;
376    JacobiRotation<Scalar> rot;
377    RealScalar n = sqrt(numext::abs2(work_matrix.coeff(p,p)) + numext::abs2(work_matrix.coeff(q,p)));
378
379    if(n==0)
380    {
381      z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
382      work_matrix.row(p) *= z;
383      if(svd.computeU()) svd.m_matrixU.col(p) *= conj(z);
384      if(work_matrix.coeff(q,q)!=Scalar(0))
385      {
386        z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
387        work_matrix.row(q) *= z;
388        if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
389      }
390      // otherwise the second row is already zero, so we have nothing to do.
391    }
392    else
393    {
394      rot.c() = conj(work_matrix.coeff(p,p)) / n;
395      rot.s() = work_matrix.coeff(q,p) / n;
396      work_matrix.applyOnTheLeft(p,q,rot);
397      if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
398      if(work_matrix.coeff(p,q) != Scalar(0))
399      {
400        Scalar z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
401        work_matrix.col(q) *= z;
402        if(svd.computeV()) svd.m_matrixV.col(q) *= z;
403      }
404      if(work_matrix.coeff(q,q) != Scalar(0))
405      {
406        z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
407        work_matrix.row(q) *= z;
408        if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
409      }
410    }
411  }
412};
413
414template<typename MatrixType, typename RealScalar, typename Index>
415void real_2x2_jacobi_svd(const MatrixType& matrix, Index p, Index q,
416                            JacobiRotation<RealScalar> *j_left,
417                            JacobiRotation<RealScalar> *j_right)
418{
419  using std::sqrt;
420  using std::abs;
421  Matrix<RealScalar,2,2> m;
422  m << numext::real(matrix.coeff(p,p)), numext::real(matrix.coeff(p,q)),
423       numext::real(matrix.coeff(q,p)), numext::real(matrix.coeff(q,q));
424  JacobiRotation<RealScalar> rot1;
425  RealScalar t = m.coeff(0,0) + m.coeff(1,1);
426  RealScalar d = m.coeff(1,0) - m.coeff(0,1);
427  if(t == RealScalar(0))
428  {
429    rot1.c() = RealScalar(0);
430    rot1.s() = d > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
431  }
432  else
433  {
434    RealScalar t2d2 = numext::hypot(t,d);
435    rot1.c() = abs(t)/t2d2;
436    rot1.s() = d/t2d2;
437    if(t<RealScalar(0))
438      rot1.s() = -rot1.s();
439  }
440  m.applyOnTheLeft(0,1,rot1);
441  j_right->makeJacobi(m,0,1);
442  *j_left  = rot1 * j_right->transpose();
443}
444
445} // end namespace internal
446
447/** \ingroup SVD_Module
448  *
449  *
450  * \class JacobiSVD
451  *
452  * \brief Two-sided Jacobi SVD decomposition of a rectangular matrix
453  *
454  * \param MatrixType the type of the matrix of which we are computing the SVD decomposition
455  * \param QRPreconditioner this optional parameter allows to specify the type of QR decomposition that will be used internally
456  *                        for the R-SVD step for non-square matrices. See discussion of possible values below.
457  *
458  * SVD decomposition consists in decomposing any n-by-p matrix \a A as a product
459  *   \f[ A = U S V^* \f]
460  * 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;
461  * 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
462  * and right \em singular \em vectors of \a A respectively.
463  *
464  * Singular values are always sorted in decreasing order.
465  *
466  * This JacobiSVD decomposition computes only the singular values by default. If you want \a U or \a V, you need to ask for them explicitly.
467  *
468  * 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
469  * 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
470  * 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,
471  * 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.
472  *
473  * Here's an example demonstrating basic usage:
474  * \include JacobiSVD_basic.cpp
475  * Output: \verbinclude JacobiSVD_basic.out
476  *
477  * This JacobiSVD class is a two-sided Jacobi R-SVD decomposition, ensuring optimal reliability and accuracy. The downside is that it's slower than
478  * bidiagonalizing SVD algorithms for large square matrices; however its complexity is still \f$ O(n^2p) \f$ where \a n is the smaller dimension and
479  * \a p is the greater dimension, meaning that it is still of the same order of complexity as the faster bidiagonalizing R-SVD algorithms.
480  * In particular, like any R-SVD, it takes advantage of non-squareness in that its complexity is only linear in the greater dimension.
481  *
482  * If the input matrix has inf or nan coefficients, the result of the computation is undefined, but the computation is guaranteed to
483  * terminate in finite (and reasonable) time.
484  *
485  * The possible values for QRPreconditioner are:
486  * \li ColPivHouseholderQRPreconditioner is the default. In practice it's very safe. It uses column-pivoting QR.
487  * \li FullPivHouseholderQRPreconditioner, is the safest and slowest. It uses full-pivoting QR.
488  *     Contrary to other QRs, it doesn't allow computing thin unitaries.
489  * \li HouseholderQRPreconditioner is the fastest, and less safe and accurate than the pivoting variants. It uses non-pivoting QR.
490  *     This is very similar in safety and accuracy to the bidiagonalization process used by bidiagonalizing SVD algorithms (since bidiagonalization
491  *     is inherently non-pivoting). However the resulting SVD is still more reliable than bidiagonalizing SVDs because the Jacobi-based iterarive
492  *     process is more reliable than the optimized bidiagonal SVD iterations.
493  * \li NoQRPreconditioner allows not to use a QR preconditioner at all. This is useful if you know that you will only be computing
494  *     JacobiSVD decompositions of square matrices. Non-square matrices require a QR preconditioner. Using this option will result in
495  *     faster compilation and smaller executable code. It won't significantly speed up computation, since JacobiSVD is always checking
496  *     if QR preconditioning is needed before applying it anyway.
497  *
498  * \sa MatrixBase::jacobiSvd()
499  */
500template<typename _MatrixType, int QRPreconditioner> class JacobiSVD
501{
502  public:
503
504    typedef _MatrixType MatrixType;
505    typedef typename MatrixType::Scalar Scalar;
506    typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
507    typedef typename MatrixType::Index Index;
508    enum {
509      RowsAtCompileTime = MatrixType::RowsAtCompileTime,
510      ColsAtCompileTime = MatrixType::ColsAtCompileTime,
511      DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
512      MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
513      MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
514      MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
515      MatrixOptions = MatrixType::Options
516    };
517
518    typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
519                   MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
520            MatrixUType;
521    typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime,
522                   MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime>
523            MatrixVType;
524    typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
525    typedef typename internal::plain_row_type<MatrixType>::type RowType;
526    typedef typename internal::plain_col_type<MatrixType>::type ColType;
527    typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
528                   MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
529            WorkMatrixType;
530
531    /** \brief Default Constructor.
532      *
533      * The default constructor is useful in cases in which the user intends to
534      * perform decompositions via JacobiSVD::compute(const MatrixType&).
535      */
536    JacobiSVD()
537      : m_isInitialized(false),
538        m_isAllocated(false),
539        m_usePrescribedThreshold(false),
540        m_computationOptions(0),
541        m_rows(-1), m_cols(-1), m_diagSize(0)
542    {}
543
544
545    /** \brief Default Constructor with memory preallocation
546      *
547      * Like the default constructor but with preallocation of the internal data
548      * according to the specified problem size.
549      * \sa JacobiSVD()
550      */
551    JacobiSVD(Index rows, Index cols, unsigned int computationOptions = 0)
552      : m_isInitialized(false),
553        m_isAllocated(false),
554        m_usePrescribedThreshold(false),
555        m_computationOptions(0),
556        m_rows(-1), m_cols(-1)
557    {
558      allocate(rows, cols, computationOptions);
559    }
560
561    /** \brief Constructor performing the decomposition of given matrix.
562     *
563     * \param matrix the matrix to decompose
564     * \param computationOptions optional parameter allowing to specify if you want full or thin U or V unitaries to be computed.
565     *                           By default, none is computed. This is a bit-field, the possible bits are #ComputeFullU, #ComputeThinU,
566     *                           #ComputeFullV, #ComputeThinV.
567     *
568     * Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not
569     * available with the (non-default) FullPivHouseholderQR preconditioner.
570     */
571    JacobiSVD(const MatrixType& matrix, unsigned int computationOptions = 0)
572      : m_isInitialized(false),
573        m_isAllocated(false),
574        m_usePrescribedThreshold(false),
575        m_computationOptions(0),
576        m_rows(-1), m_cols(-1)
577    {
578      compute(matrix, computationOptions);
579    }
580
581    /** \brief Method performing the decomposition of given matrix using custom options.
582     *
583     * \param matrix the matrix to decompose
584     * \param computationOptions optional parameter allowing to specify if you want full or thin U or V unitaries to be computed.
585     *                           By default, none is computed. This is a bit-field, the possible bits are #ComputeFullU, #ComputeThinU,
586     *                           #ComputeFullV, #ComputeThinV.
587     *
588     * Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not
589     * available with the (non-default) FullPivHouseholderQR preconditioner.
590     */
591    JacobiSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
592
593    /** \brief Method performing the decomposition of given matrix using current options.
594     *
595     * \param matrix the matrix to decompose
596     *
597     * This method uses the current \a computationOptions, as already passed to the constructor or to compute(const MatrixType&, unsigned int).
598     */
599    JacobiSVD& compute(const MatrixType& matrix)
600    {
601      return compute(matrix, m_computationOptions);
602    }
603
604    /** \returns the \a U matrix.
605     *
606     * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p,
607     * the U matrix is n-by-n if you asked for #ComputeFullU, and is n-by-m if you asked for #ComputeThinU.
608     *
609     * The \a m first columns of \a U are the left singular vectors of the matrix being decomposed.
610     *
611     * This method asserts that you asked for \a U to be computed.
612     */
613    const MatrixUType& matrixU() const
614    {
615      eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
616      eigen_assert(computeU() && "This JacobiSVD decomposition didn't compute U. Did you ask for it?");
617      return m_matrixU;
618    }
619
620    /** \returns the \a V matrix.
621     *
622     * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p,
623     * the V matrix is p-by-p if you asked for #ComputeFullV, and is p-by-m if you asked for ComputeThinV.
624     *
625     * The \a m first columns of \a V are the right singular vectors of the matrix being decomposed.
626     *
627     * This method asserts that you asked for \a V to be computed.
628     */
629    const MatrixVType& matrixV() const
630    {
631      eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
632      eigen_assert(computeV() && "This JacobiSVD decomposition didn't compute V. Did you ask for it?");
633      return m_matrixV;
634    }
635
636    /** \returns the vector of singular values.
637     *
638     * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p, the
639     * returned vector has size \a m.  Singular values are always sorted in decreasing order.
640     */
641    const SingularValuesType& singularValues() const
642    {
643      eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
644      return m_singularValues;
645    }
646
647    /** \returns true if \a U (full or thin) is asked for in this SVD decomposition */
648    inline bool computeU() const { return m_computeFullU || m_computeThinU; }
649    /** \returns true if \a V (full or thin) is asked for in this SVD decomposition */
650    inline bool computeV() const { return m_computeFullV || m_computeThinV; }
651
652    /** \returns a (least squares) solution of \f$ A x = b \f$ using the current SVD decomposition of A.
653      *
654      * \param b the right-hand-side of the equation to solve.
655      *
656      * \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.
657      *
658      * \note SVD solving is implicitly least-squares. Thus, this method serves both purposes of exact solving and least-squares solving.
659      * In other words, the returned solution is guaranteed to minimize the Euclidean norm \f$ \Vert A x - b \Vert \f$.
660      */
661    template<typename Rhs>
662    inline const internal::solve_retval<JacobiSVD, Rhs>
663    solve(const MatrixBase<Rhs>& b) const
664    {
665      eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
666      eigen_assert(computeU() && computeV() && "JacobiSVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
667      return internal::solve_retval<JacobiSVD, Rhs>(*this, b.derived());
668    }
669
670    /** \returns the number of singular values that are not exactly 0 */
671    Index nonzeroSingularValues() const
672    {
673      eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
674      return m_nonzeroSingularValues;
675    }
676
677    /** \returns the rank of the matrix of which \c *this is the SVD.
678      *
679      * \note This method has to determine which singular values should be considered nonzero.
680      *       For that, it uses the threshold value that you can control by calling
681      *       setThreshold(const RealScalar&).
682      */
683    inline Index rank() const
684    {
685      using std::abs;
686      eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
687      if(m_singularValues.size()==0) return 0;
688      RealScalar premultiplied_threshold = m_singularValues.coeff(0) * threshold();
689      Index i = m_nonzeroSingularValues-1;
690      while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
691      return i+1;
692    }
693
694    /** Allows to prescribe a threshold to be used by certain methods, such as rank() and solve(),
695      * which need to determine when singular values are to be considered nonzero.
696      * This is not used for the SVD decomposition itself.
697      *
698      * When it needs to get the threshold value, Eigen calls threshold().
699      * The default is \c NumTraits<Scalar>::epsilon()
700      *
701      * \param threshold The new value to use as the threshold.
702      *
703      * A singular value will be considered nonzero if its value is strictly greater than
704      *  \f$ \vert singular value \vert \leqslant threshold \times \vert max singular value \vert \f$.
705      *
706      * If you want to come back to the default behavior, call setThreshold(Default_t)
707      */
708    JacobiSVD& setThreshold(const RealScalar& threshold)
709    {
710      m_usePrescribedThreshold = true;
711      m_prescribedThreshold = threshold;
712      return *this;
713    }
714
715    /** Allows to come back to the default behavior, letting Eigen use its default formula for
716      * determining the threshold.
717      *
718      * You should pass the special object Eigen::Default as parameter here.
719      * \code svd.setThreshold(Eigen::Default); \endcode
720      *
721      * See the documentation of setThreshold(const RealScalar&).
722      */
723    JacobiSVD& setThreshold(Default_t)
724    {
725      m_usePrescribedThreshold = false;
726      return *this;
727    }
728
729    /** Returns the threshold that will be used by certain methods such as rank().
730      *
731      * See the documentation of setThreshold(const RealScalar&).
732      */
733    RealScalar threshold() const
734    {
735      eigen_assert(m_isInitialized || m_usePrescribedThreshold);
736      return m_usePrescribedThreshold ? m_prescribedThreshold
737                                      : (std::max<Index>)(1,m_diagSize)*NumTraits<Scalar>::epsilon();
738    }
739
740    inline Index rows() const { return m_rows; }
741    inline Index cols() const { return m_cols; }
742
743  private:
744    void allocate(Index rows, Index cols, unsigned int computationOptions);
745
746  protected:
747    MatrixUType m_matrixU;
748    MatrixVType m_matrixV;
749    SingularValuesType m_singularValues;
750    WorkMatrixType m_workMatrix;
751    bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
752    bool m_computeFullU, m_computeThinU;
753    bool m_computeFullV, m_computeThinV;
754    unsigned int m_computationOptions;
755    Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
756    RealScalar m_prescribedThreshold;
757
758    template<typename __MatrixType, int _QRPreconditioner, bool _IsComplex>
759    friend struct internal::svd_precondition_2x2_block_to_be_real;
760    template<typename __MatrixType, int _QRPreconditioner, int _Case, bool _DoAnything>
761    friend struct internal::qr_preconditioner_impl;
762
763    internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols;
764    internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows;
765};
766
767template<typename MatrixType, int QRPreconditioner>
768void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols, unsigned int computationOptions)
769{
770  eigen_assert(rows >= 0 && cols >= 0);
771
772  if (m_isAllocated &&
773      rows == m_rows &&
774      cols == m_cols &&
775      computationOptions == m_computationOptions)
776  {
777    return;
778  }
779
780  m_rows = rows;
781  m_cols = cols;
782  m_isInitialized = false;
783  m_isAllocated = true;
784  m_computationOptions = computationOptions;
785  m_computeFullU = (computationOptions & ComputeFullU) != 0;
786  m_computeThinU = (computationOptions & ComputeThinU) != 0;
787  m_computeFullV = (computationOptions & ComputeFullV) != 0;
788  m_computeThinV = (computationOptions & ComputeThinV) != 0;
789  eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U");
790  eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V");
791  eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
792              "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
793  if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
794  {
795      eigen_assert(!(m_computeThinU || m_computeThinV) &&
796              "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
797              "Use the ColPivHouseholderQR preconditioner instead.");
798  }
799  m_diagSize = (std::min)(m_rows, m_cols);
800  m_singularValues.resize(m_diagSize);
801  if(RowsAtCompileTime==Dynamic)
802    m_matrixU.resize(m_rows, m_computeFullU ? m_rows
803                            : m_computeThinU ? m_diagSize
804                            : 0);
805  if(ColsAtCompileTime==Dynamic)
806    m_matrixV.resize(m_cols, m_computeFullV ? m_cols
807                            : m_computeThinV ? m_diagSize
808                            : 0);
809  m_workMatrix.resize(m_diagSize, m_diagSize);
810
811  if(m_cols>m_rows) m_qr_precond_morecols.allocate(*this);
812  if(m_rows>m_cols) m_qr_precond_morerows.allocate(*this);
813}
814
815template<typename MatrixType, int QRPreconditioner>
816JacobiSVD<MatrixType, QRPreconditioner>&
817JacobiSVD<MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)
818{
819  using std::abs;
820  allocate(matrix.rows(), matrix.cols(), computationOptions);
821
822  // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
823  // only worsening the precision of U and V as we accumulate more rotations
824  const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
825
826  // limit for very small denormal numbers to be considered zero in order to avoid infinite loops (see bug 286)
827  const RealScalar considerAsZero = RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min();
828
829  /*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */
830
831  if(!m_qr_precond_morecols.run(*this, matrix) && !m_qr_precond_morerows.run(*this, matrix))
832  {
833    m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize);
834    if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
835    if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
836    if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
837    if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
838  }
839
840  // Scaling factor to reduce over/under-flows
841  RealScalar scale = m_workMatrix.cwiseAbs().maxCoeff();
842  if(scale==RealScalar(0)) scale = RealScalar(1);
843  m_workMatrix /= scale;
844
845  /*** step 2. The main Jacobi SVD iteration. ***/
846
847  bool finished = false;
848  while(!finished)
849  {
850    finished = true;
851
852    // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
853
854    for(Index p = 1; p < m_diagSize; ++p)
855    {
856      for(Index q = 0; q < p; ++q)
857      {
858        // if this 2x2 sub-matrix is not diagonal already...
859        // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
860        // keep us iterating forever. Similarly, small denormal numbers are considered zero.
861        using std::max;
862        RealScalar threshold = (max)(considerAsZero, precision * (max)(abs(m_workMatrix.coeff(p,p)),
863                                                                       abs(m_workMatrix.coeff(q,q))));
864        if((max)(abs(m_workMatrix.coeff(p,q)),abs(m_workMatrix.coeff(q,p))) > threshold)
865        {
866          finished = false;
867
868          // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
869          internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q);
870          JacobiRotation<RealScalar> j_left, j_right;
871          internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
872
873          // accumulate resulting Jacobi rotations
874          m_workMatrix.applyOnTheLeft(p,q,j_left);
875          if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
876
877          m_workMatrix.applyOnTheRight(p,q,j_right);
878          if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
879        }
880      }
881    }
882  }
883
884  /*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/
885
886  for(Index i = 0; i < m_diagSize; ++i)
887  {
888    RealScalar a = abs(m_workMatrix.coeff(i,i));
889    m_singularValues.coeffRef(i) = a;
890    if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
891  }
892
893  /*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/
894
895  m_nonzeroSingularValues = m_diagSize;
896  for(Index i = 0; i < m_diagSize; i++)
897  {
898    Index pos;
899    RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
900    if(maxRemainingSingularValue == RealScalar(0))
901    {
902      m_nonzeroSingularValues = i;
903      break;
904    }
905    if(pos)
906    {
907      pos += i;
908      std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
909      if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
910      if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
911    }
912  }
913
914  m_singularValues *= scale;
915
916  m_isInitialized = true;
917  return *this;
918}
919
920namespace internal {
921template<typename _MatrixType, int QRPreconditioner, typename Rhs>
922struct solve_retval<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
923  : solve_retval_base<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
924{
925  typedef JacobiSVD<_MatrixType, QRPreconditioner> JacobiSVDType;
926  EIGEN_MAKE_SOLVE_HELPERS(JacobiSVDType,Rhs)
927
928  template<typename Dest> void evalTo(Dest& dst) const
929  {
930    eigen_assert(rhs().rows() == dec().rows());
931
932    // A = U S V^*
933    // So A^{-1} = V S^{-1} U^*
934
935    Matrix<Scalar, Dynamic, Rhs::ColsAtCompileTime, 0, _MatrixType::MaxRowsAtCompileTime, Rhs::MaxColsAtCompileTime> tmp;
936    Index rank = dec().rank();
937
938    tmp.noalias() = dec().matrixU().leftCols(rank).adjoint() * rhs();
939    tmp = dec().singularValues().head(rank).asDiagonal().inverse() * tmp;
940    dst = dec().matrixV().leftCols(rank) * tmp;
941  }
942};
943} // end namespace internal
944
945/** \svd_module
946  *
947  * \return the singular value decomposition of \c *this computed by two-sided
948  * Jacobi transformations.
949  *
950  * \sa class JacobiSVD
951  */
952template<typename Derived>
953JacobiSVD<typename MatrixBase<Derived>::PlainObject>
954MatrixBase<Derived>::jacobiSvd(unsigned int computationOptions) const
955{
956  return JacobiSVD<PlainObject>(*this, computationOptions);
957}
958
959} // end namespace Eigen
960
961#endif // EIGEN_JACOBISVD_H
962