Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
JacobiSVD.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009-2010 Benoit Jacob <[email protected]>
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 
13 namespace Eigen {
14 
15 namespace internal {
16 // forward declaration (needed by ICC)
17 // the empty body is required by MSVC
18 template<typename MatrixType, int QRPreconditioner,
19  bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
20 struct 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 
29 enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
30 
31 template<typename MatrixType, int QRPreconditioner, int Case>
32 struct 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 
46 template<typename MatrixType, int QRPreconditioner, int Case,
47  bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret
48 > struct qr_preconditioner_impl {};
49 
50 template<typename MatrixType, int QRPreconditioner, int Case>
51 class qr_preconditioner_impl<MatrixType, QRPreconditioner, Case, false>
52 {
53 public:
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 
64 template<typename MatrixType>
65 class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
66 {
67 public:
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  }
99 private:
100  typedef FullPivHouseholderQR<MatrixType> QRType;
101  QRType m_qr;
102  WorkspaceType m_workspace;
103 };
104 
105 template<typename MatrixType>
106 class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
107 {
108 public:
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  }
146 private:
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 
155 template<typename MatrixType>
156 class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
157 {
158 public:
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 
190 private:
191  typedef ColPivHouseholderQR<MatrixType> QRType;
192  QRType m_qr;
193  typename internal::plain_col_type<MatrixType>::type m_workspace;
194 };
195 
196 template<typename MatrixType>
197 class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
198 {
199 public:
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 
246 private:
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 
255 template<typename MatrixType>
256 class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
257 {
258 public:
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  }
289 private:
290  typedef HouseholderQR<MatrixType> QRType;
291  QRType m_qr;
292  typename internal::plain_col_type<MatrixType>::type m_workspace;
293 };
294 
295 template<typename MatrixType>
296 class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
297 {
298 public:
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 
345 private:
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 
357 template<typename MatrixType, int QRPreconditioner>
358 struct 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 
365 template<typename MatrixType, int QRPreconditioner>
366 struct 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 
414 template<typename MatrixType, typename RealScalar, typename Index>
415 void 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 
500 template<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 
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 
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 
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 
591  JacobiSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
592 
599  JacobiSVD& compute(const MatrixType& matrix)
600  {
601  return compute(matrix, m_computationOptions);
602  }
603 
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 
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 
641  const SingularValuesType& singularValues() const
642  {
643  eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
644  return m_singularValues;
645  }
646 
648  inline bool computeU() const { return m_computeFullU || m_computeThinU; }
650  inline bool computeV() const { return m_computeFullV || m_computeThinV; }
651 
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 
671  Index nonzeroSingularValues() const
672  {
673  eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
674  return m_nonzeroSingularValues;
675  }
676 
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 
708  JacobiSVD& setThreshold(const RealScalar& threshold)
709  {
710  m_usePrescribedThreshold = true;
711  m_prescribedThreshold = threshold;
712  return *this;
713  }
714 
724  {
725  m_usePrescribedThreshold = false;
726  return *this;
727  }
728 
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  static void check_template_parameters()
747  {
748  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
749  }
750 
751  protected:
752  MatrixUType m_matrixU;
753  MatrixVType m_matrixV;
754  SingularValuesType m_singularValues;
755  WorkMatrixType m_workMatrix;
756  bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
757  bool m_computeFullU, m_computeThinU;
758  bool m_computeFullV, m_computeThinV;
759  unsigned int m_computationOptions;
760  Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
761  RealScalar m_prescribedThreshold;
762 
763  template<typename __MatrixType, int _QRPreconditioner, bool _IsComplex>
764  friend struct internal::svd_precondition_2x2_block_to_be_real;
765  template<typename __MatrixType, int _QRPreconditioner, int _Case, bool _DoAnything>
766  friend struct internal::qr_preconditioner_impl;
767 
768  internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols;
769  internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows;
770  MatrixType m_scaledMatrix;
771 };
772 
773 template<typename MatrixType, int QRPreconditioner>
774 void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols, unsigned int computationOptions)
775 {
776  eigen_assert(rows >= 0 && cols >= 0);
777 
778  if (m_isAllocated &&
779  rows == m_rows &&
780  cols == m_cols &&
781  computationOptions == m_computationOptions)
782  {
783  return;
784  }
785 
786  m_rows = rows;
787  m_cols = cols;
788  m_isInitialized = false;
789  m_isAllocated = true;
790  m_computationOptions = computationOptions;
791  m_computeFullU = (computationOptions & ComputeFullU) != 0;
792  m_computeThinU = (computationOptions & ComputeThinU) != 0;
793  m_computeFullV = (computationOptions & ComputeFullV) != 0;
794  m_computeThinV = (computationOptions & ComputeThinV) != 0;
795  eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U");
796  eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V");
797  eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
798  "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
799  if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
800  {
801  eigen_assert(!(m_computeThinU || m_computeThinV) &&
802  "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
803  "Use the ColPivHouseholderQR preconditioner instead.");
804  }
805  m_diagSize = (std::min)(m_rows, m_cols);
806  m_singularValues.resize(m_diagSize);
807  if(RowsAtCompileTime==Dynamic)
808  m_matrixU.resize(m_rows, m_computeFullU ? m_rows
809  : m_computeThinU ? m_diagSize
810  : 0);
811  if(ColsAtCompileTime==Dynamic)
812  m_matrixV.resize(m_cols, m_computeFullV ? m_cols
813  : m_computeThinV ? m_diagSize
814  : 0);
815  m_workMatrix.resize(m_diagSize, m_diagSize);
816 
817  if(m_cols>m_rows) m_qr_precond_morecols.allocate(*this);
818  if(m_rows>m_cols) m_qr_precond_morerows.allocate(*this);
819  if(m_rows!=m_cols) m_scaledMatrix.resize(rows,cols);
820 }
821 
822 template<typename MatrixType, int QRPreconditioner>
823 JacobiSVD<MatrixType, QRPreconditioner>&
824 JacobiSVD<MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)
825 {
826  check_template_parameters();
827 
828  using std::abs;
829  allocate(matrix.rows(), matrix.cols(), computationOptions);
830 
831  // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
832  // only worsening the precision of U and V as we accumulate more rotations
833  const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
834 
835  // limit for very small denormal numbers to be considered zero in order to avoid infinite loops (see bug 286)
836  const RealScalar considerAsZero = RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min();
837 
838  // Scaling factor to reduce over/under-flows
839  RealScalar scale = matrix.cwiseAbs().maxCoeff();
840  if(scale==RealScalar(0)) scale = RealScalar(1);
841 
842  /*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */
843 
844  if(m_rows!=m_cols)
845  {
846  m_scaledMatrix = matrix / scale;
847  m_qr_precond_morecols.run(*this, m_scaledMatrix);
848  m_qr_precond_morerows.run(*this, m_scaledMatrix);
849  }
850  else
851  {
852  m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize) / scale;
853  if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
854  if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
855  if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
856  if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
857  }
858 
859  /*** step 2. The main Jacobi SVD iteration. ***/
860 
861  bool finished = false;
862  while(!finished)
863  {
864  finished = true;
865 
866  // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
867 
868  for(Index p = 1; p < m_diagSize; ++p)
869  {
870  for(Index q = 0; q < p; ++q)
871  {
872  // if this 2x2 sub-matrix is not diagonal already...
873  // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
874  // keep us iterating forever. Similarly, small denormal numbers are considered zero.
875  using std::max;
876  RealScalar threshold = (max)(considerAsZero, precision * (max)(abs(m_workMatrix.coeff(p,p)),
877  abs(m_workMatrix.coeff(q,q))));
878  // We compare both values to threshold instead of calling max to be robust to NaN (See bug 791)
879  if(abs(m_workMatrix.coeff(p,q))>threshold || abs(m_workMatrix.coeff(q,p)) > threshold)
880  {
881  finished = false;
882 
883  // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
884  internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q);
885  JacobiRotation<RealScalar> j_left, j_right;
886  internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
887 
888  // accumulate resulting Jacobi rotations
889  m_workMatrix.applyOnTheLeft(p,q,j_left);
890  if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
891 
892  m_workMatrix.applyOnTheRight(p,q,j_right);
893  if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
894  }
895  }
896  }
897  }
898 
899  /*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/
900 
901  for(Index i = 0; i < m_diagSize; ++i)
902  {
903  RealScalar a = abs(m_workMatrix.coeff(i,i));
904  m_singularValues.coeffRef(i) = a;
905  if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
906  }
907 
908  /*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/
909 
910  m_nonzeroSingularValues = m_diagSize;
911  for(Index i = 0; i < m_diagSize; i++)
912  {
913  Index pos;
914  RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
915  if(maxRemainingSingularValue == RealScalar(0))
916  {
917  m_nonzeroSingularValues = i;
918  break;
919  }
920  if(pos)
921  {
922  pos += i;
923  std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
924  if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
925  if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
926  }
927  }
928 
929  m_singularValues *= scale;
930 
931  m_isInitialized = true;
932  return *this;
933 }
934 
935 namespace internal {
936 template<typename _MatrixType, int QRPreconditioner, typename Rhs>
937 struct solve_retval<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
938  : solve_retval_base<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
939 {
940  typedef JacobiSVD<_MatrixType, QRPreconditioner> JacobiSVDType;
941  EIGEN_MAKE_SOLVE_HELPERS(JacobiSVDType,Rhs)
942 
943  template<typename Dest> void evalTo(Dest& dst) const
944  {
945  eigen_assert(rhs().rows() == dec().rows());
946 
947  // A = U S V^*
948  // So A^{-1} = V S^{-1} U^*
949 
951  Index rank = dec().rank();
952 
953  tmp.noalias() = dec().matrixU().leftCols(rank).adjoint() * rhs();
954  tmp = dec().singularValues().head(rank).asDiagonal().inverse() * tmp;
955  dst = dec().matrixV().leftCols(rank) * tmp;
956  }
957 };
958 } // end namespace internal
959 
967 template<typename Derived>
968 JacobiSVD<typename MatrixBase<Derived>::PlainObject>
969 MatrixBase<Derived>::jacobiSvd(unsigned int computationOptions) const
970 {
971  return JacobiSVD<PlainObject>(*this, computationOptions);
972 }
973 
974 } // end namespace Eigen
975 
976 #endif // EIGEN_JACOBISVD_H
Index rank() const
Definition: JacobiSVD.h:683
JacobiSVD< PlainObject > jacobiSvd(unsigned int computationOptions=0) const
Definition: JacobiSVD.h:969
Definition: Constants.h:359
bool computeV() const
Definition: JacobiSVD.h:650
Definition: Constants.h:329
const internal::solve_retval< JacobiSVD, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: JacobiSVD.h:663
RealScalar threshold() const
Definition: JacobiSVD.h:733
NoAlias< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols >, Eigen::MatrixBase > noalias()
Rotation given by a cosine-sine pair.
Definition: ForwardDeclarations.h:228
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const int Dynamic
Definition: Constants.h:21
Definition: Constants.h:331
const SingularValuesType & singularValues() const
Definition: JacobiSVD.h:641
ColsBlockXpr leftCols(Index n)
Definition: DenseBase.h:527
SegmentReturnType head(Index n)
Definition: DenseBase.h:806
JacobiRotation transpose() const
Definition: Jacobi.h:59
JacobiSVD()
Default Constructor.
Definition: JacobiSVD.h:536
const MatrixUType & matrixU() const
Definition: JacobiSVD.h:613
JacobiSVD(Index rows, Index cols, unsigned int computationOptions=0)
Default Constructor with memory preallocation.
Definition: JacobiSVD.h:551
JacobiSVD & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
Definition: JacobiSVD.h:824
Definition: Constants.h:333
Index nonzeroSingularValues() const
Definition: JacobiSVD.h:671
const MatrixVType & matrixV() const
Definition: JacobiSVD.h:629
bool computeU() const
Definition: JacobiSVD.h:648
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:224
JacobiSVD & setThreshold(Default_t)
Definition: JacobiSVD.h:723
JacobiSVD & setThreshold(const RealScalar &threshold)
Definition: JacobiSVD.h:708
JacobiSVD(const MatrixType &matrix, unsigned int computationOptions=0)
Constructor performing the decomposition of given matrix.
Definition: JacobiSVD.h:571
Definition: Constants.h:327
Definition: Constants.h:361
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Definition: Constants.h:363
JacobiSVD & compute(const MatrixType &matrix)
Method performing the decomposition of given matrix using current options.
Definition: JacobiSVD.h:599