$darkmode
Eigen  5.0.1-dev
BDCSVD.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // We used the "A Divide-And-Conquer Algorithm for the Bidiagonal SVD"
5 // research report written by Ming Gu and Stanley C.Eisenstat
6 // The code variable names correspond to the names they used in their
7 // report
8 //
9 // Copyright (C) 2013 Gauthier Brun <brun.gauthier@gmail.com>
10 // Copyright (C) 2013 Nicolas Carre <nicolas.carre@ensimag.fr>
11 // Copyright (C) 2013 Jean Ceccato <jean.ceccato@ensimag.fr>
12 // Copyright (C) 2013 Pierre Zoppitelli <pierre.zoppitelli@ensimag.fr>
13 // Copyright (C) 2013 Jitse Niesen <jitse@maths.leeds.ac.uk>
14 // Copyright (C) 2014-2017 Gael Guennebaud <gael.guennebaud@inria.fr>
15 //
16 // Source Code Form is subject to the terms of the Mozilla
17 // Public License v. 2.0. If a copy of the MPL was not distributed
18 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
19 
20 #ifndef EIGEN_BDCSVD_H
21 #define EIGEN_BDCSVD_H
22 // #define EIGEN_BDCSVD_DEBUG_VERBOSE
23 // #define EIGEN_BDCSVD_SANITY_CHECKS
24 
25 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
26 #undef eigen_internal_assert
27 #define eigen_internal_assert(X) assert(X);
28 #endif
29 
30 // IWYU pragma: private
31 #include "./InternalHeaderCheck.h"
32 
33 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
34 #include <iostream>
35 #endif
36 
37 namespace Eigen {
38 
39 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
40 IOFormat bdcsvdfmt(8, 0, ", ", "\n", " [", "]");
41 #endif
42 
43 template <typename MatrixType_, int Options>
44 class BDCSVD;
45 
46 namespace internal {
47 
48 template <typename MatrixType_, int Options>
49 struct traits<BDCSVD<MatrixType_, Options> > : svd_traits<MatrixType_, Options> {
50  typedef MatrixType_ MatrixType;
51 };
52 
53 } // end namespace internal
54 
84 template <typename MatrixType_, int Options_>
85 class BDCSVD : public SVDBase<BDCSVD<MatrixType_, Options_> > {
86  typedef SVDBase<BDCSVD> Base;
87 
88  public:
89  using Base::cols;
90  using Base::computeU;
91  using Base::computeV;
92  using Base::diagSize;
93  using Base::rows;
94 
95  typedef MatrixType_ MatrixType;
96  typedef typename Base::Scalar Scalar;
97  typedef typename Base::RealScalar RealScalar;
98  typedef typename NumTraits<RealScalar>::Literal Literal;
99  typedef typename Base::Index Index;
100  enum {
101  Options = Options_,
102  QRDecomposition = Options & internal::QRPreconditionerBits,
103  ComputationOptions = Options & internal::ComputationOptionsBits,
104  RowsAtCompileTime = Base::RowsAtCompileTime,
105  ColsAtCompileTime = Base::ColsAtCompileTime,
106  DiagSizeAtCompileTime = Base::DiagSizeAtCompileTime,
107  MaxRowsAtCompileTime = Base::MaxRowsAtCompileTime,
108  MaxColsAtCompileTime = Base::MaxColsAtCompileTime,
109  MaxDiagSizeAtCompileTime = Base::MaxDiagSizeAtCompileTime,
110  MatrixOptions = Base::MatrixOptions
111  };
112 
113  typedef typename Base::MatrixUType MatrixUType;
114  typedef typename Base::MatrixVType MatrixVType;
115  typedef typename Base::SingularValuesType SingularValuesType;
116 
117  typedef Matrix<Scalar, Dynamic, Dynamic, ColMajor> MatrixX;
118  typedef Matrix<RealScalar, Dynamic, Dynamic, ColMajor> MatrixXr;
119  typedef Matrix<RealScalar, Dynamic, 1> VectorType;
120  typedef Array<RealScalar, Dynamic, 1> ArrayXr;
121  typedef Array<Index, 1, Dynamic> ArrayXi;
122  typedef Ref<ArrayXr> ArrayRef;
123  typedef Ref<ArrayXi> IndicesRef;
124 
130  BDCSVD() : m_algoswap(16), m_isTranspose(false), m_compU(false), m_compV(false), m_numIters(0) {}
131 
138  BDCSVD(Index rows, Index cols) : m_algoswap(16), m_numIters(0) {
139  allocate(rows, cols, internal::get_computation_options(Options));
140  }
141 
158  EIGEN_DEPRECATED_WITH_REASON("Options should be specified using the class template parameter.")
159  BDCSVD(Index rows, Index cols, unsigned int computationOptions) : m_algoswap(16), m_numIters(0) {
160  internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, rows, cols);
161  allocate(rows, cols, computationOptions);
162  }
163 
169  template <typename Derived>
170  BDCSVD(const MatrixBase<Derived>& matrix) : m_algoswap(16), m_numIters(0) {
171  compute_impl(matrix, internal::get_computation_options(Options));
172  }
173 
186  template <typename Derived>
187  EIGEN_DEPRECATED_WITH_REASON("Options should be specified using the class template parameter.")
188  BDCSVD(const MatrixBase<Derived>& matrix, unsigned int computationOptions) : m_algoswap(16), m_numIters(0) {
189  internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, matrix.rows(), matrix.cols());
190  compute_impl(matrix, computationOptions);
191  }
192 
193  ~BDCSVD() {}
194 
200  template <typename Derived>
202  return compute_impl(matrix, m_computationOptions);
203  }
204 
214  template <typename Derived>
215  EIGEN_DEPRECATED_WITH_REASON("Options should be specified using the class template parameter.")
216  BDCSVD& compute(const MatrixBase<Derived>& matrix, unsigned int computationOptions) {
217  internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, matrix.rows(), matrix.cols());
218  return compute_impl(matrix, computationOptions);
219  }
220 
221  void setSwitchSize(int s) {
222  eigen_assert(s >= 3 && "BDCSVD the size of the algo switch has to be at least 3.");
223  m_algoswap = s;
224  }
225 
226  private:
227  template <typename Derived>
228  BDCSVD& compute_impl(const MatrixBase<Derived>& matrix, unsigned int computationOptions);
229  void divide(Index firstCol, Index lastCol, Index firstRowW, Index firstColW, Index shift);
230  void computeSVDofM(Index firstCol, Index n, MatrixXr& U, VectorType& singVals, MatrixXr& V);
231  void computeSingVals(const ArrayRef& col0, const ArrayRef& diag, const IndicesRef& perm, VectorType& singVals,
232  ArrayRef shifts, ArrayRef mus);
233  void perturbCol0(const ArrayRef& col0, const ArrayRef& diag, const IndicesRef& perm, const VectorType& singVals,
234  const ArrayRef& shifts, const ArrayRef& mus, ArrayRef zhat);
235  void computeSingVecs(const ArrayRef& zhat, const ArrayRef& diag, const IndicesRef& perm, const VectorType& singVals,
236  const ArrayRef& shifts, const ArrayRef& mus, MatrixXr& U, MatrixXr& V);
237  void deflation43(Index firstCol, Index shift, Index i, Index size);
238  void deflation44(Index firstColu, Index firstColm, Index firstRowW, Index firstColW, Index i, Index j, Index size);
239  void deflation(Index firstCol, Index lastCol, Index k, Index firstRowW, Index firstColW, Index shift);
240  template <typename HouseholderU, typename HouseholderV, typename NaiveU, typename NaiveV>
241  void copyUV(const HouseholderU& householderU, const HouseholderV& householderV, const NaiveU& naiveU,
242  const NaiveV& naivev);
243  void structured_update(Block<MatrixXr, Dynamic, Dynamic> A, const MatrixXr& B, Index n1);
244  static RealScalar secularEq(RealScalar x, const ArrayRef& col0, const ArrayRef& diag, const IndicesRef& perm,
245  const ArrayRef& diagShifted, RealScalar shift);
246  template <typename SVDType>
247  void computeBaseCase(SVDType& svd, Index n, Index firstCol, Index firstRowW, Index firstColW, Index shift);
248 
249  protected:
250  void allocate(Index rows, Index cols, unsigned int computationOptions);
251  MatrixXr m_naiveU, m_naiveV;
252  MatrixXr m_computed;
253  Index m_nRec;
254  ArrayXr m_workspace;
255  ArrayXi m_workspaceI;
256  int m_algoswap;
257  bool m_isTranspose, m_compU, m_compV, m_useQrDecomp;
258  JacobiSVD<MatrixType, ComputationOptions> smallSvd;
259  HouseholderQR<MatrixX> qrDecomp;
260  internal::UpperBidiagonalization<MatrixX> bid;
261  MatrixX copyWorkspace;
262  MatrixX reducedTriangle;
263 
264  using Base::m_computationOptions;
265  using Base::m_computeThinU;
266  using Base::m_computeThinV;
267  using Base::m_info;
268  using Base::m_isInitialized;
269  using Base::m_matrixU;
270  using Base::m_matrixV;
271  using Base::m_nonzeroSingularValues;
272  using Base::m_singularValues;
273 
274  public:
275  int m_numIters;
276 }; // end class BDCSVD
277 
278 // Method to allocate and initialize matrix and attributes
279 template <typename MatrixType, int Options>
280 void BDCSVD<MatrixType, Options>::allocate(Index rows, Index cols, unsigned int computationOptions) {
281  if (Base::allocate(rows, cols, computationOptions)) return;
282 
283  if (cols < m_algoswap)
284  smallSvd.allocate(rows, cols, Options == 0 ? computationOptions : internal::get_computation_options(Options));
285 
286  m_computed = MatrixXr::Zero(diagSize() + 1, diagSize());
287  m_compU = computeV();
288  m_compV = computeU();
289  m_isTranspose = (cols > rows);
290  if (m_isTranspose) std::swap(m_compU, m_compV);
291 
292  // kMinAspectRatio is the crossover point that determines if we perform R-Bidiagonalization
293  // or bidiagonalize the input matrix directly.
294  // It is based off of LAPACK's dgesdd routine, which uses 11.0/6.0
295  // we use a larger scalar to prevent a regression for relatively square matrices.
296  constexpr Index kMinAspectRatio = 4;
297  constexpr bool disableQrDecomp = static_cast<int>(QRDecomposition) == static_cast<int>(DisableQRDecomposition);
298  m_useQrDecomp = !disableQrDecomp && ((rows / kMinAspectRatio > cols) || (cols / kMinAspectRatio > rows));
299  if (m_useQrDecomp) {
300  qrDecomp = HouseholderQR<MatrixX>((std::max)(rows, cols), (std::min)(rows, cols));
301  reducedTriangle = MatrixX(diagSize(), diagSize());
302  }
303 
304  copyWorkspace = MatrixX(m_isTranspose ? cols : rows, m_isTranspose ? rows : cols);
305  bid = internal::UpperBidiagonalization<MatrixX>(m_useQrDecomp ? diagSize() : copyWorkspace.rows(),
306  m_useQrDecomp ? diagSize() : copyWorkspace.cols());
307 
308  if (m_compU)
309  m_naiveU = MatrixXr::Zero(diagSize() + 1, diagSize() + 1);
310  else
311  m_naiveU = MatrixXr::Zero(2, diagSize() + 1);
312 
313  if (m_compV) m_naiveV = MatrixXr::Zero(diagSize(), diagSize());
314 
315  m_workspace.resize((diagSize() + 1) * (diagSize() + 1) * 3);
316  m_workspaceI.resize(3 * diagSize());
317 } // end allocate
318 
319 template <typename MatrixType, int Options>
320 template <typename Derived>
321 BDCSVD<MatrixType, Options>& BDCSVD<MatrixType, Options>::compute_impl(const MatrixBase<Derived>& matrix,
322  unsigned int computationOptions) {
323  EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived, MatrixType);
324  EIGEN_STATIC_ASSERT((std::is_same<typename Derived::Scalar, typename MatrixType::Scalar>::value),
325  Input matrix must have the same Scalar type as the BDCSVD object.);
326 
327 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
328  std::cout << "\n\n\n================================================================================================="
329  "=====================\n\n\n";
330 #endif
331  using std::abs;
332 
333  allocate(matrix.rows(), matrix.cols(), computationOptions);
334 
335  const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
336 
337  //**** step -1 - If the problem is too small, directly falls back to JacobiSVD and return
338  if (matrix.cols() < m_algoswap) {
339  smallSvd.compute(matrix);
340  m_isInitialized = true;
341  m_info = smallSvd.info();
342  if (m_info == Success || m_info == NoConvergence) {
343  if (computeU()) m_matrixU = smallSvd.matrixU();
344  if (computeV()) m_matrixV = smallSvd.matrixV();
345  m_singularValues = smallSvd.singularValues();
346  m_nonzeroSingularValues = smallSvd.nonzeroSingularValues();
347  }
348  return *this;
349  }
350 
351  //**** step 0 - Copy the input matrix and apply scaling to reduce over/under-flows
352  RealScalar scale = matrix.cwiseAbs().template maxCoeff<PropagateNaN>();
353  if (!(numext::isfinite)(scale)) {
354  m_isInitialized = true;
355  m_info = InvalidInput;
356  return *this;
357  }
358 
359  if (numext::is_exactly_zero(scale)) scale = Literal(1);
360 
361  if (m_isTranspose)
362  copyWorkspace = matrix.adjoint() / scale;
363  else
364  copyWorkspace = matrix / scale;
365 
366  //**** step 1 - Bidiagonalization.
367  // If the problem is sufficiently rectangular, we perform R-Bidiagonalization: compute A = Q(R/0)
368  // and then bidiagonalize R. Otherwise, if the problem is relatively square, we
369  // bidiagonalize the input matrix directly.
370  if (m_useQrDecomp) {
371  qrDecomp.compute(copyWorkspace);
372  reducedTriangle = qrDecomp.matrixQR().topRows(diagSize());
373  reducedTriangle.template triangularView<StrictlyLower>().setZero();
374  bid.compute(reducedTriangle);
375  } else {
376  bid.compute(copyWorkspace);
377  }
378 
379  //**** step 2 - Divide & Conquer
380  m_naiveU.setZero();
381  m_naiveV.setZero();
382  // FIXME this line involves a temporary matrix
383  m_computed.topRows(diagSize()) = bid.bidiagonal().toDenseMatrix().transpose();
384  m_computed.template bottomRows<1>().setZero();
385  divide(0, diagSize() - 1, 0, 0, 0);
386  if (m_info != Success && m_info != NoConvergence) {
387  m_isInitialized = true;
388  return *this;
389  }
390 
391  //**** step 3 - Copy singular values and vectors
392  for (int i = 0; i < diagSize(); i++) {
393  RealScalar a = abs(m_computed.coeff(i, i));
394  m_singularValues.coeffRef(i) = a * scale;
395  if (a < considerZero) {
396  m_nonzeroSingularValues = i;
397  m_singularValues.tail(diagSize() - i - 1).setZero();
398  break;
399  } else if (i == diagSize() - 1) {
400  m_nonzeroSingularValues = i + 1;
401  break;
402  }
403  }
404 
405  //**** step 4 - Finalize unitaries U and V
406  if (m_isTranspose)
407  copyUV(bid.householderV(), bid.householderU(), m_naiveV, m_naiveU);
408  else
409  copyUV(bid.householderU(), bid.householderV(), m_naiveU, m_naiveV);
410 
411  if (m_useQrDecomp) {
412  if (m_isTranspose && computeV())
413  m_matrixV.applyOnTheLeft(qrDecomp.householderQ());
414  else if (!m_isTranspose && computeU())
415  m_matrixU.applyOnTheLeft(qrDecomp.householderQ());
416  }
417 
418  m_isInitialized = true;
419  return *this;
420 } // end compute
421 
422 template <typename MatrixType, int Options>
423 template <typename HouseholderU, typename HouseholderV, typename NaiveU, typename NaiveV>
424 void BDCSVD<MatrixType, Options>::copyUV(const HouseholderU& householderU, const HouseholderV& householderV,
425  const NaiveU& naiveU, const NaiveV& naiveV) {
426  // Note exchange of U and V: m_matrixU is set from m_naiveV and vice versa
427  if (computeU()) {
428  Index Ucols = m_computeThinU ? diagSize() : rows();
429  m_matrixU = MatrixX::Identity(rows(), Ucols);
430  m_matrixU.topLeftCorner(diagSize(), diagSize()) =
431  naiveV.template cast<Scalar>().topLeftCorner(diagSize(), diagSize());
432  // FIXME the following conditionals involve temporary buffers
433  if (m_useQrDecomp)
434  m_matrixU.topLeftCorner(householderU.cols(), diagSize()).applyOnTheLeft(householderU);
435  else
436  m_matrixU.applyOnTheLeft(householderU);
437  }
438  if (computeV()) {
439  Index Vcols = m_computeThinV ? diagSize() : cols();
440  m_matrixV = MatrixX::Identity(cols(), Vcols);
441  m_matrixV.topLeftCorner(diagSize(), diagSize()) =
442  naiveU.template cast<Scalar>().topLeftCorner(diagSize(), diagSize());
443  // FIXME the following conditionals involve temporary buffers
444  if (m_useQrDecomp)
445  m_matrixV.topLeftCorner(householderV.cols(), diagSize()).applyOnTheLeft(householderV);
446  else
447  m_matrixV.applyOnTheLeft(householderV);
448  }
449 }
450 
459 template <typename MatrixType, int Options>
460 void BDCSVD<MatrixType, Options>::structured_update(Block<MatrixXr, Dynamic, Dynamic> A, const MatrixXr& B, Index n1) {
461  Index n = A.rows();
462  if (n > 100) {
463  // If the matrices are large enough, let's exploit the sparse structure of A by
464  // splitting it in half (wrt n1), and packing the non-zero columns.
465  Index n2 = n - n1;
466  Map<MatrixXr> A1(m_workspace.data(), n1, n);
467  Map<MatrixXr> A2(m_workspace.data() + n1 * n, n2, n);
468  Map<MatrixXr> B1(m_workspace.data() + n * n, n, n);
469  Map<MatrixXr> B2(m_workspace.data() + 2 * n * n, n, n);
470  Index k1 = 0, k2 = 0;
471  for (Index j = 0; j < n; ++j) {
472  if ((A.col(j).head(n1).array() != Literal(0)).any()) {
473  A1.col(k1) = A.col(j).head(n1);
474  B1.row(k1) = B.row(j);
475  ++k1;
476  }
477  if ((A.col(j).tail(n2).array() != Literal(0)).any()) {
478  A2.col(k2) = A.col(j).tail(n2);
479  B2.row(k2) = B.row(j);
480  ++k2;
481  }
482  }
483 
484  A.topRows(n1).noalias() = A1.leftCols(k1) * B1.topRows(k1);
485  A.bottomRows(n2).noalias() = A2.leftCols(k2) * B2.topRows(k2);
486  } else {
487  Map<MatrixXr, Aligned> tmp(m_workspace.data(), n, n);
488  tmp.noalias() = A * B;
489  A = tmp;
490  }
491 }
492 
493 template <typename MatrixType, int Options>
494 template <typename SVDType>
495 void BDCSVD<MatrixType, Options>::computeBaseCase(SVDType& svd, Index n, Index firstCol, Index firstRowW,
496  Index firstColW, Index shift) {
497  svd.compute(m_computed.block(firstCol, firstCol, n + 1, n));
498  m_info = svd.info();
499  if (m_info != Success && m_info != NoConvergence) return;
500  if (m_compU)
501  m_naiveU.block(firstCol, firstCol, n + 1, n + 1).real() = svd.matrixU();
502  else {
503  m_naiveU.row(0).segment(firstCol, n + 1).real() = svd.matrixU().row(0);
504  m_naiveU.row(1).segment(firstCol, n + 1).real() = svd.matrixU().row(n);
505  }
506  if (m_compV) m_naiveV.block(firstRowW, firstColW, n, n).real() = svd.matrixV();
507  m_computed.block(firstCol + shift, firstCol + shift, n + 1, n).setZero();
508  m_computed.diagonal().segment(firstCol + shift, n) = svd.singularValues().head(n);
509 }
510 
511 // The divide algorithm is done "in place", we are always working on subsets of the same matrix. The divide methods
512 // takes as argument the place of the submatrix we are currently working on.
513 
514 //@param firstCol : The Index of the first column of the submatrix of m_computed and for m_naiveU;
515 //@param lastCol : The Index of the last column of the submatrix of m_computed and for m_naiveU;
516 // lastCol + 1 - firstCol is the size of the submatrix.
517 //@param firstRowW : The Index of the first row of the matrix W that we are to change. (see the reference paper section
518 // 1 for more information on W)
519 //@param firstColW : Same as firstRowW with the column.
520 //@param shift : Each time one takes the left submatrix, one must add 1 to the shift. Why? Because! We actually want the
521 // last column of the U submatrix
522 // to become the first column (*coeff) and to shift all the other columns to the right. There are more details on the
523 // reference paper.
524 template <typename MatrixType, int Options>
525 void BDCSVD<MatrixType, Options>::divide(Index firstCol, Index lastCol, Index firstRowW, Index firstColW, Index shift) {
526  // requires rows = cols + 1;
527  using std::abs;
528  using std::pow;
529  using std::sqrt;
530  const Index n = lastCol - firstCol + 1;
531  const Index k = n / 2;
532  const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
533  RealScalar alphaK;
534  RealScalar betaK;
535  RealScalar r0;
536  RealScalar lambda, phi, c0, s0;
537  VectorType l, f;
538  // We use the other algorithm which is more efficient for small
539  // matrices.
540  if (n < m_algoswap) {
541  // FIXME this block involves temporaries
542  if (m_compV) {
543  JacobiSVD<MatrixXr, ComputeFullU | ComputeFullV> baseSvd;
544  computeBaseCase(baseSvd, n, firstCol, firstRowW, firstColW, shift);
545  } else {
546  JacobiSVD<MatrixXr, ComputeFullU> baseSvd;
547  computeBaseCase(baseSvd, n, firstCol, firstRowW, firstColW, shift);
548  }
549  return;
550  }
551  // We use the divide and conquer algorithm
552  alphaK = m_computed(firstCol + k, firstCol + k);
553  betaK = m_computed(firstCol + k + 1, firstCol + k);
554  // The divide must be done in that order in order to have good results. Divide change the data inside the submatrices
555  // and the divide of the right submatrice reads one column of the left submatrice. That's why we need to treat the
556  // right submatrix before the left one.
557  divide(k + 1 + firstCol, lastCol, k + 1 + firstRowW, k + 1 + firstColW, shift);
558  if (m_info != Success && m_info != NoConvergence) return;
559  divide(firstCol, k - 1 + firstCol, firstRowW, firstColW + 1, shift + 1);
560  if (m_info != Success && m_info != NoConvergence) return;
561 
562  if (m_compU) {
563  lambda = m_naiveU(firstCol + k, firstCol + k);
564  phi = m_naiveU(firstCol + k + 1, lastCol + 1);
565  } else {
566  lambda = m_naiveU(1, firstCol + k);
567  phi = m_naiveU(0, lastCol + 1);
568  }
569  r0 = sqrt((abs(alphaK * lambda) * abs(alphaK * lambda)) + abs(betaK * phi) * abs(betaK * phi));
570  if (m_compU) {
571  l = m_naiveU.row(firstCol + k).segment(firstCol, k);
572  f = m_naiveU.row(firstCol + k + 1).segment(firstCol + k + 1, n - k - 1);
573  } else {
574  l = m_naiveU.row(1).segment(firstCol, k);
575  f = m_naiveU.row(0).segment(firstCol + k + 1, n - k - 1);
576  }
577  if (m_compV) m_naiveV(firstRowW + k, firstColW) = Literal(1);
578  if (r0 < considerZero) {
579  c0 = Literal(1);
580  s0 = Literal(0);
581  } else {
582  c0 = alphaK * lambda / r0;
583  s0 = betaK * phi / r0;
584  }
585 
586 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
587  eigen_internal_assert(m_naiveU.allFinite());
588  eigen_internal_assert(m_naiveV.allFinite());
589  eigen_internal_assert(m_computed.allFinite());
590 #endif
591 
592  if (m_compU) {
593  MatrixXr q1(m_naiveU.col(firstCol + k).segment(firstCol, k + 1));
594  // we shiftW Q1 to the right
595  for (Index i = firstCol + k - 1; i >= firstCol; i--)
596  m_naiveU.col(i + 1).segment(firstCol, k + 1) = m_naiveU.col(i).segment(firstCol, k + 1);
597  // we shift q1 at the left with a factor c0
598  m_naiveU.col(firstCol).segment(firstCol, k + 1) = (q1 * c0);
599  // last column = q1 * - s0
600  m_naiveU.col(lastCol + 1).segment(firstCol, k + 1) = (q1 * (-s0));
601  // first column = q2 * s0
602  m_naiveU.col(firstCol).segment(firstCol + k + 1, n - k) =
603  m_naiveU.col(lastCol + 1).segment(firstCol + k + 1, n - k) * s0;
604  // q2 *= c0
605  m_naiveU.col(lastCol + 1).segment(firstCol + k + 1, n - k) *= c0;
606  } else {
607  RealScalar q1 = m_naiveU(0, firstCol + k);
608  // we shift Q1 to the right
609  for (Index i = firstCol + k - 1; i >= firstCol; i--) m_naiveU(0, i + 1) = m_naiveU(0, i);
610  // we shift q1 at the left with a factor c0
611  m_naiveU(0, firstCol) = (q1 * c0);
612  // last column = q1 * - s0
613  m_naiveU(0, lastCol + 1) = (q1 * (-s0));
614  // first column = q2 * s0
615  m_naiveU(1, firstCol) = m_naiveU(1, lastCol + 1) * s0;
616  // q2 *= c0
617  m_naiveU(1, lastCol + 1) *= c0;
618  m_naiveU.row(1).segment(firstCol + 1, k).setZero();
619  m_naiveU.row(0).segment(firstCol + k + 1, n - k - 1).setZero();
620  }
621 
622 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
623  eigen_internal_assert(m_naiveU.allFinite());
624  eigen_internal_assert(m_naiveV.allFinite());
625  eigen_internal_assert(m_computed.allFinite());
626 #endif
627 
628  m_computed(firstCol + shift, firstCol + shift) = r0;
629  m_computed.col(firstCol + shift).segment(firstCol + shift + 1, k) = alphaK * l.transpose().real();
630  m_computed.col(firstCol + shift).segment(firstCol + shift + k + 1, n - k - 1) = betaK * f.transpose().real();
631 
632 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
633  ArrayXr tmp1 = (m_computed.block(firstCol + shift, firstCol + shift, n, n)).jacobiSvd().singularValues();
634 #endif
635  // Second part: try to deflate singular values in combined matrix
636  deflation(firstCol, lastCol, k, firstRowW, firstColW, shift);
637 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
638  ArrayXr tmp2 = (m_computed.block(firstCol + shift, firstCol + shift, n, n)).jacobiSvd().singularValues();
639  std::cout << "\n\nj1 = " << tmp1.transpose().format(bdcsvdfmt) << "\n";
640  std::cout << "j2 = " << tmp2.transpose().format(bdcsvdfmt) << "\n\n";
641  std::cout << "err: " << ((tmp1 - tmp2).abs() > 1e-12 * tmp2.abs()).transpose() << "\n";
642  static int count = 0;
643  std::cout << "# " << ++count << "\n\n";
644  eigen_internal_assert((tmp1 - tmp2).matrix().norm() < 1e-14 * tmp2.matrix().norm());
645 // eigen_internal_assert(count<681);
646 // eigen_internal_assert(((tmp1-tmp2).abs()<1e-13*tmp2.abs()).all());
647 #endif
648 
649  // Third part: compute SVD of combined matrix
650  MatrixXr UofSVD, VofSVD;
651  VectorType singVals;
652  computeSVDofM(firstCol + shift, n, UofSVD, singVals, VofSVD);
653 
654 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
655  eigen_internal_assert(UofSVD.allFinite());
656  eigen_internal_assert(VofSVD.allFinite());
657 #endif
658 
659  if (m_compU)
660  structured_update(m_naiveU.block(firstCol, firstCol, n + 1, n + 1), UofSVD, (n + 2) / 2);
661  else {
662  Map<Matrix<RealScalar, 2, Dynamic>, Aligned> tmp(m_workspace.data(), 2, n + 1);
663  tmp.noalias() = m_naiveU.middleCols(firstCol, n + 1) * UofSVD;
664  m_naiveU.middleCols(firstCol, n + 1) = tmp;
665  }
666 
667  if (m_compV) structured_update(m_naiveV.block(firstRowW, firstColW, n, n), VofSVD, (n + 1) / 2);
668 
669 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
670  eigen_internal_assert(m_naiveU.allFinite());
671  eigen_internal_assert(m_naiveV.allFinite());
672  eigen_internal_assert(m_computed.allFinite());
673 #endif
674 
675  m_computed.block(firstCol + shift, firstCol + shift, n, n).setZero();
676  m_computed.block(firstCol + shift, firstCol + shift, n, n).diagonal() = singVals;
677 } // end divide
678 
679 // Compute SVD of m_computed.block(firstCol, firstCol, n + 1, n); this block only has non-zeros in
680 // the first column and on the diagonal and has undergone deflation, so diagonal is in increasing
681 // order except for possibly the (0,0) entry. The computed SVD is stored U, singVals and V, except
682 // that if m_compV is false, then V is not computed. Singular values are sorted in decreasing order.
683 //
684 // TODO Opportunities for optimization: better root finding algo, better stopping criterion, better
685 // handling of round-off errors, be consistent in ordering
686 // For instance, to solve the secular equation using FMM, see
687 // http://www.stat.uchicago.edu/~lekheng/courses/302/classics/greengard-rokhlin.pdf
688 template <typename MatrixType, int Options>
689 void BDCSVD<MatrixType, Options>::computeSVDofM(Index firstCol, Index n, MatrixXr& U, VectorType& singVals,
690  MatrixXr& V) {
691  const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
692  using std::abs;
693  ArrayRef col0 = m_computed.col(firstCol).segment(firstCol, n);
694  m_workspace.head(n) = m_computed.block(firstCol, firstCol, n, n).diagonal();
695  ArrayRef diag = m_workspace.head(n);
696  diag(0) = Literal(0);
697 
698  // Allocate space for singular values and vectors
699  singVals.resize(n);
700  U.resize(n + 1, n + 1);
701  if (m_compV) V.resize(n, n);
702 
703 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
704  if (col0.hasNaN() || diag.hasNaN()) std::cout << "\n\nHAS NAN\n\n";
705 #endif
706 
707  // Many singular values might have been deflated, the zero ones have been moved to the end,
708  // but others are interleaved and we must ignore them at this stage.
709  // To this end, let's compute a permutation skipping them:
710  Index actual_n = n;
711  while (actual_n > 1 && numext::is_exactly_zero(diag(actual_n - 1))) {
712  --actual_n;
713  eigen_internal_assert(numext::is_exactly_zero(col0(actual_n)));
714  }
715  Index m = 0; // size of the deflated problem
716  for (Index k = 0; k < actual_n; ++k)
717  if (abs(col0(k)) > considerZero) m_workspaceI(m++) = k;
718  Map<ArrayXi> perm(m_workspaceI.data(), m);
719 
720  Map<ArrayXr> shifts(m_workspace.data() + 1 * n, n);
721  Map<ArrayXr> mus(m_workspace.data() + 2 * n, n);
722  Map<ArrayXr> zhat(m_workspace.data() + 3 * n, n);
723 
724 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
725  std::cout << "computeSVDofM using:\n";
726  std::cout << " z: " << col0.transpose() << "\n";
727  std::cout << " d: " << diag.transpose() << "\n";
728 #endif
729 
730  // Compute singVals, shifts, and mus
731  computeSingVals(col0, diag, perm, singVals, shifts, mus);
732 
733 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
734  std::cout << " j: "
735  << (m_computed.block(firstCol, firstCol, n, n)).jacobiSvd().singularValues().transpose().reverse()
736  << "\n\n";
737  std::cout << " sing-val: " << singVals.transpose() << "\n";
738  std::cout << " mu: " << mus.transpose() << "\n";
739  std::cout << " shift: " << shifts.transpose() << "\n";
740 
741  {
742  std::cout << "\n\n mus: " << mus.head(actual_n).transpose() << "\n\n";
743  std::cout << " check1 (expect0) : "
744  << ((singVals.array() - (shifts + mus)) / singVals.array()).head(actual_n).transpose() << "\n\n";
745  eigen_internal_assert((((singVals.array() - (shifts + mus)) / singVals.array()).head(actual_n) >= 0).all());
746  std::cout << " check2 (>0) : " << ((singVals.array() - diag) / singVals.array()).head(actual_n).transpose()
747  << "\n\n";
748  eigen_internal_assert((((singVals.array() - diag) / singVals.array()).head(actual_n) >= 0).all());
749  }
750 #endif
751 
752 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
753  eigen_internal_assert(singVals.allFinite());
754  eigen_internal_assert(mus.allFinite());
755  eigen_internal_assert(shifts.allFinite());
756 #endif
757 
758  // Compute zhat
759  perturbCol0(col0, diag, perm, singVals, shifts, mus, zhat);
760 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
761  std::cout << " zhat: " << zhat.transpose() << "\n";
762 #endif
763 
764 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
765  eigen_internal_assert(zhat.allFinite());
766 #endif
767 
768  computeSingVecs(zhat, diag, perm, singVals, shifts, mus, U, V);
769 
770 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
771  std::cout << "U^T U: " << (U.transpose() * U - MatrixXr(MatrixXr::Identity(U.cols(), U.cols()))).norm() << "\n";
772  std::cout << "V^T V: " << (V.transpose() * V - MatrixXr(MatrixXr::Identity(V.cols(), V.cols()))).norm() << "\n";
773 #endif
774 
775 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
776  eigen_internal_assert(m_naiveU.allFinite());
777  eigen_internal_assert(m_naiveV.allFinite());
778  eigen_internal_assert(m_computed.allFinite());
779  eigen_internal_assert(U.allFinite());
780  eigen_internal_assert(V.allFinite());
781 // eigen_internal_assert((U.transpose() * U - MatrixXr(MatrixXr::Identity(U.cols(),U.cols()))).norm() <
782 // 100*NumTraits<RealScalar>::epsilon() * n); eigen_internal_assert((V.transpose() * V -
783 // MatrixXr(MatrixXr::Identity(V.cols(),V.cols()))).norm() < 100*NumTraits<RealScalar>::epsilon() * n);
784 #endif
785 
786  // Because of deflation, the singular values might not be completely sorted.
787  // Fortunately, reordering them is a O(n) problem
788  for (Index i = 0; i < actual_n - 1; ++i) {
789  if (singVals(i) > singVals(i + 1)) {
790  using std::swap;
791  swap(singVals(i), singVals(i + 1));
792  U.col(i).swap(U.col(i + 1));
793  if (m_compV) V.col(i).swap(V.col(i + 1));
794  }
795  }
796 
797 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
798  {
799  bool singular_values_sorted =
800  (((singVals.segment(1, actual_n - 1) - singVals.head(actual_n - 1))).array() >= 0).all();
801  if (!singular_values_sorted)
802  std::cout << "Singular values are not sorted: " << singVals.segment(1, actual_n).transpose() << "\n";
803  eigen_internal_assert(singular_values_sorted);
804  }
805 #endif
806 
807  // Reverse order so that singular values in increased order
808  // Because of deflation, the zeros singular-values are already at the end
809  singVals.head(actual_n).reverseInPlace();
810  U.leftCols(actual_n).rowwise().reverseInPlace();
811  if (m_compV) V.leftCols(actual_n).rowwise().reverseInPlace();
812 
813 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
814  JacobiSVD<MatrixXr> jsvd(m_computed.block(firstCol, firstCol, n, n));
815  std::cout << " * j: " << jsvd.singularValues().transpose() << "\n\n";
816  std::cout << " * sing-val: " << singVals.transpose() << "\n";
817 // std::cout << " * err: " << ((jsvd.singularValues()-singVals)>1e-13*singVals.norm()).transpose() << "\n";
818 #endif
819 }
820 
821 template <typename MatrixType, int Options>
822 typename BDCSVD<MatrixType, Options>::RealScalar BDCSVD<MatrixType, Options>::secularEq(
823  RealScalar mu, const ArrayRef& col0, const ArrayRef& diag, const IndicesRef& perm, const ArrayRef& diagShifted,
824  RealScalar shift) {
825  Index m = perm.size();
826  RealScalar res = Literal(1);
827  for (Index i = 0; i < m; ++i) {
828  Index j = perm(i);
829  // The following expression could be rewritten to involve only a single division,
830  // but this would make the expression more sensitive to overflow.
831  res += (col0(j) / (diagShifted(j) - mu)) * (col0(j) / (diag(j) + shift + mu));
832  }
833  return res;
834 }
835 
836 template <typename MatrixType, int Options>
837 void BDCSVD<MatrixType, Options>::computeSingVals(const ArrayRef& col0, const ArrayRef& diag, const IndicesRef& perm,
838  VectorType& singVals, ArrayRef shifts, ArrayRef mus) {
839  using std::abs;
840  using std::sqrt;
841  using std::swap;
842 
843  Index n = col0.size();
844  Index actual_n = n;
845  // Note that here actual_n is computed based on col0(i)==0 instead of diag(i)==0 as above
846  // because 1) we have diag(i)==0 => col0(i)==0 and 2) if col0(i)==0, then diag(i) is already a singular value.
847  while (actual_n > 1 && numext::is_exactly_zero(col0(actual_n - 1))) --actual_n;
848 
849  for (Index k = 0; k < n; ++k) {
850  if (numext::is_exactly_zero(col0(k)) || actual_n == 1) {
851  // if col0(k) == 0, then entry is deflated, so singular value is on diagonal
852  // if actual_n==1, then the deflated problem is already diagonalized
853  singVals(k) = k == 0 ? col0(0) : diag(k);
854  mus(k) = Literal(0);
855  shifts(k) = k == 0 ? col0(0) : diag(k);
856  continue;
857  }
858 
859  // otherwise, use secular equation to find singular value
860  RealScalar left = diag(k);
861  RealScalar right; // was: = (k != actual_n-1) ? diag(k+1) : (diag(actual_n-1) + col0.matrix().norm());
862  if (k == actual_n - 1)
863  right = (diag(actual_n - 1) + col0.matrix().norm());
864  else {
865  // Skip deflated singular values,
866  // recall that at this stage we assume that z[j]!=0 and all entries for which z[j]==0 have been put aside.
867  // This should be equivalent to using perm[]
868  Index l = k + 1;
869  while (numext::is_exactly_zero(col0(l))) {
870  ++l;
871  eigen_internal_assert(l < actual_n);
872  }
873  right = diag(l);
874  }
875 
876  // first decide whether it's closer to the left end or the right end
877  RealScalar mid = left + (right - left) / Literal(2);
878  RealScalar fMid = secularEq(mid, col0, diag, perm, diag, Literal(0));
879 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
880  std::cout << "right-left = " << right - left << "\n";
881  // std::cout << "fMid = " << fMid << " " << secularEq(mid-left, col0, diag, perm, ArrayXr(diag-left), left)
882  // << " " << secularEq(mid-right, col0, diag, perm, ArrayXr(diag-right), right) <<
883  // "\n";
884  std::cout << " = " << secularEq(left + RealScalar(0.000001) * (right - left), col0, diag, perm, diag, 0) << " "
885  << secularEq(left + RealScalar(0.1) * (right - left), col0, diag, perm, diag, 0) << " "
886  << secularEq(left + RealScalar(0.2) * (right - left), col0, diag, perm, diag, 0) << " "
887  << secularEq(left + RealScalar(0.3) * (right - left), col0, diag, perm, diag, 0) << " "
888  << secularEq(left + RealScalar(0.4) * (right - left), col0, diag, perm, diag, 0) << " "
889  << secularEq(left + RealScalar(0.49) * (right - left), col0, diag, perm, diag, 0) << " "
890  << secularEq(left + RealScalar(0.5) * (right - left), col0, diag, perm, diag, 0) << " "
891  << secularEq(left + RealScalar(0.51) * (right - left), col0, diag, perm, diag, 0) << " "
892  << secularEq(left + RealScalar(0.6) * (right - left), col0, diag, perm, diag, 0) << " "
893  << secularEq(left + RealScalar(0.7) * (right - left), col0, diag, perm, diag, 0) << " "
894  << secularEq(left + RealScalar(0.8) * (right - left), col0, diag, perm, diag, 0) << " "
895  << secularEq(left + RealScalar(0.9) * (right - left), col0, diag, perm, diag, 0) << " "
896  << secularEq(left + RealScalar(0.999999) * (right - left), col0, diag, perm, diag, 0) << "\n";
897 #endif
898  RealScalar shift = (k == actual_n - 1 || fMid > Literal(0)) ? left : right;
899 
900  // measure everything relative to shift
901  Map<ArrayXr> diagShifted(m_workspace.data() + 4 * n, n);
902  diagShifted = diag - shift;
903 
904  if (k != actual_n - 1) {
905  // check that after the shift, f(mid) is still negative:
906  RealScalar midShifted = (right - left) / RealScalar(2);
907  // we can test exact equality here, because shift comes from `... ? left : right`
908  if (numext::equal_strict(shift, right)) midShifted = -midShifted;
909  RealScalar fMidShifted = secularEq(midShifted, col0, diag, perm, diagShifted, shift);
910  if (fMidShifted > 0) {
911  // fMid was erroneous, fix it:
912  shift = fMidShifted > Literal(0) ? left : right;
913  diagShifted = diag - shift;
914  }
915  }
916 
917  // initial guess
918  RealScalar muPrev, muCur;
919  // we can test exact equality here, because shift comes from `... ? left : right`
920  if (numext::equal_strict(shift, left)) {
921  muPrev = (right - left) * RealScalar(0.1);
922  if (k == actual_n - 1)
923  muCur = right - left;
924  else
925  muCur = (right - left) * RealScalar(0.5);
926  } else {
927  muPrev = -(right - left) * RealScalar(0.1);
928  muCur = -(right - left) * RealScalar(0.5);
929  }
930 
931  RealScalar fPrev = secularEq(muPrev, col0, diag, perm, diagShifted, shift);
932  RealScalar fCur = secularEq(muCur, col0, diag, perm, diagShifted, shift);
933  if (abs(fPrev) < abs(fCur)) {
934  swap(fPrev, fCur);
935  swap(muPrev, muCur);
936  }
937 
938  // rational interpolation: fit a function of the form a / mu + b through the two previous
939  // iterates and use its zero to compute the next iterate
940  bool useBisection = fPrev * fCur > Literal(0);
941  while (!numext::is_exactly_zero(fCur) &&
942  abs(muCur - muPrev) >
943  Literal(8) * NumTraits<RealScalar>::epsilon() * numext::maxi<RealScalar>(abs(muCur), abs(muPrev)) &&
944  abs(fCur - fPrev) > NumTraits<RealScalar>::epsilon() && !useBisection) {
945  ++m_numIters;
946 
947  // Find a and b such that the function f(mu) = a / mu + b matches the current and previous samples.
948  RealScalar a = (fCur - fPrev) / (Literal(1) / muCur - Literal(1) / muPrev);
949  RealScalar b = fCur - a / muCur;
950  // And find mu such that f(mu)==0:
951  RealScalar muZero = -a / b;
952  RealScalar fZero = secularEq(muZero, col0, diag, perm, diagShifted, shift);
953 
954 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
955  eigen_internal_assert((numext::isfinite)(fZero));
956 #endif
957 
958  muPrev = muCur;
959  fPrev = fCur;
960  muCur = muZero;
961  fCur = fZero;
962 
963  // we can test exact equality here, because shift comes from `... ? left : right`
964  if (numext::equal_strict(shift, left) && (muCur < Literal(0) || muCur > right - left)) useBisection = true;
965  if (numext::equal_strict(shift, right) && (muCur < -(right - left) || muCur > Literal(0))) useBisection = true;
966  if (abs(fCur) > abs(fPrev)) useBisection = true;
967  }
968 
969  // fall back on bisection method if rational interpolation did not work
970  if (useBisection) {
971 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
972  std::cout << "useBisection for k = " << k << ", actual_n = " << actual_n << "\n";
973 #endif
974  RealScalar leftShifted, rightShifted;
975  // we can test exact equality here, because shift comes from `... ? left : right`
976  if (numext::equal_strict(shift, left)) {
977  // to avoid overflow, we must have mu > max(real_min, |z(k)|/sqrt(real_max)),
978  // the factor 2 is to be more conservative
979  leftShifted =
980  numext::maxi<RealScalar>((std::numeric_limits<RealScalar>::min)(),
981  Literal(2) * abs(col0(k)) / sqrt((std::numeric_limits<RealScalar>::max)()));
982 
983  // check that we did it right:
984  eigen_internal_assert(
985  (numext::isfinite)((col0(k) / leftShifted) * (col0(k) / (diag(k) + shift + leftShifted))));
986  // I don't understand why the case k==0 would be special there:
987  // if (k == 0) rightShifted = right - left; else
988  rightShifted = (k == actual_n - 1)
989  ? right
990  : ((right - left) * RealScalar(0.51)); // theoretically we can take 0.5, but let's be safe
991  } else {
992  leftShifted = -(right - left) * RealScalar(0.51);
993  if (k + 1 < n)
994  rightShifted = -numext::maxi<RealScalar>((std::numeric_limits<RealScalar>::min)(),
995  abs(col0(k + 1)) / sqrt((std::numeric_limits<RealScalar>::max)()));
996  else
997  rightShifted = -(std::numeric_limits<RealScalar>::min)();
998  }
999 
1000  RealScalar fLeft = secularEq(leftShifted, col0, diag, perm, diagShifted, shift);
1001  eigen_internal_assert(fLeft < Literal(0));
1002 
1003 #if defined EIGEN_BDCSVD_DEBUG_VERBOSE || defined EIGEN_BDCSVD_SANITY_CHECKS || defined EIGEN_INTERNAL_DEBUGGING
1004  RealScalar fRight = secularEq(rightShifted, col0, diag, perm, diagShifted, shift);
1005 #endif
1006 
1007 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1008  if (!(numext::isfinite)(fLeft))
1009  std::cout << "f(" << leftShifted << ") =" << fLeft << " ; " << left << " " << shift << " " << right << "\n";
1010  eigen_internal_assert((numext::isfinite)(fLeft));
1011 
1012  if (!(numext::isfinite)(fRight))
1013  std::cout << "f(" << rightShifted << ") =" << fRight << " ; " << left << " " << shift << " " << right << "\n";
1014  // eigen_internal_assert((numext::isfinite)(fRight));
1015 #endif
1016 
1017 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1018  if (!(fLeft * fRight < 0)) {
1019  std::cout << "f(leftShifted) using leftShifted=" << leftShifted
1020  << " ; diagShifted(1:10):" << diagShifted.head(10).transpose() << "\n ; "
1021  << "left==shift=" << bool(left == shift) << " ; left-shift = " << (left - shift) << "\n";
1022  std::cout << "k=" << k << ", " << fLeft << " * " << fRight << " == " << fLeft * fRight << " ; "
1023  << "[" << left << " .. " << right << "] -> [" << leftShifted << " " << rightShifted
1024  << "], shift=" << shift << " , f(right)=" << secularEq(0, col0, diag, perm, diagShifted, shift)
1025  << " == " << secularEq(right, col0, diag, perm, diag, 0) << " == " << fRight << "\n";
1026  }
1027 #endif
1028  eigen_internal_assert(fLeft * fRight < Literal(0));
1029 
1030  if (fLeft < Literal(0)) {
1031  while (rightShifted - leftShifted > Literal(2) * NumTraits<RealScalar>::epsilon() *
1032  numext::maxi<RealScalar>(abs(leftShifted), abs(rightShifted))) {
1033  RealScalar midShifted = (leftShifted + rightShifted) / Literal(2);
1034  fMid = secularEq(midShifted, col0, diag, perm, diagShifted, shift);
1035  eigen_internal_assert((numext::isfinite)(fMid));
1036 
1037  if (fLeft * fMid < Literal(0)) {
1038  rightShifted = midShifted;
1039  } else {
1040  leftShifted = midShifted;
1041  fLeft = fMid;
1042  }
1043  }
1044  muCur = (leftShifted + rightShifted) / Literal(2);
1045  } else {
1046  // We have a problem as shifting on the left or right give either a positive or negative value
1047  // at the middle of [left,right]...
1048  // Instead of abbording or entering an infinite loop,
1049  // let's just use the middle as the estimated zero-crossing:
1050  muCur = (right - left) * RealScalar(0.5);
1051  // we can test exact equality here, because shift comes from `... ? left : right`
1052  if (numext::equal_strict(shift, right)) muCur = -muCur;
1053  }
1054  }
1055 
1056  singVals[k] = shift + muCur;
1057  shifts[k] = shift;
1058  mus[k] = muCur;
1059 
1060 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1061  if (k + 1 < n)
1062  std::cout << "found " << singVals[k] << " == " << shift << " + " << muCur << " from " << diag(k) << " .. "
1063  << diag(k + 1) << "\n";
1064 #endif
1065 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1066  eigen_internal_assert(k == 0 || singVals[k] >= singVals[k - 1]);
1067  eigen_internal_assert(singVals[k] >= diag(k));
1068 #endif
1069 
1070  // perturb singular value slightly if it equals diagonal entry to avoid division by zero later
1071  // (deflation is supposed to avoid this from happening)
1072  // - this does no seem to be necessary anymore -
1073  // if (singVals[k] == left) singVals[k] *= 1 + NumTraits<RealScalar>::epsilon();
1074  // if (singVals[k] == right) singVals[k] *= 1 - NumTraits<RealScalar>::epsilon();
1075  }
1076 }
1077 
1078 // zhat is perturbation of col0 for which singular vectors can be computed stably (see Section 3.1)
1079 template <typename MatrixType, int Options>
1080 void BDCSVD<MatrixType, Options>::perturbCol0(const ArrayRef& col0, const ArrayRef& diag, const IndicesRef& perm,
1081  const VectorType& singVals, const ArrayRef& shifts, const ArrayRef& mus,
1082  ArrayRef zhat) {
1083  using std::sqrt;
1084  Index n = col0.size();
1085  Index m = perm.size();
1086  if (m == 0) {
1087  zhat.setZero();
1088  return;
1089  }
1090  Index lastIdx = perm(m - 1);
1091  // The offset permits to skip deflated entries while computing zhat
1092  for (Index k = 0; k < n; ++k) {
1093  if (numext::is_exactly_zero(col0(k))) // deflated
1094  zhat(k) = Literal(0);
1095  else {
1096  // see equation (3.6)
1097  RealScalar dk = diag(k);
1098  RealScalar prod = (singVals(lastIdx) + dk) * (mus(lastIdx) + (shifts(lastIdx) - dk));
1099 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1100  if (prod < 0) {
1101  std::cout << "k = " << k << " ; z(k)=" << col0(k) << ", diag(k)=" << dk << "\n";
1102  std::cout << "prod = "
1103  << "(" << singVals(lastIdx) << " + " << dk << ") * (" << mus(lastIdx) << " + (" << shifts(lastIdx)
1104  << " - " << dk << "))"
1105  << "\n";
1106  std::cout << " = " << singVals(lastIdx) + dk << " * " << mus(lastIdx) + (shifts(lastIdx) - dk) << "\n";
1107  }
1108  eigen_internal_assert(prod >= 0);
1109 #endif
1110 
1111  for (Index l = 0; l < m; ++l) {
1112  Index i = perm(l);
1113  if (i != k) {
1114 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1115  if (i >= k && (l == 0 || l - 1 >= m)) {
1116  std::cout << "Error in perturbCol0\n";
1117  std::cout << " " << k << "/" << n << " " << l << "/" << m << " " << i << "/" << n << " ; " << col0(k)
1118  << " " << diag(k) << " "
1119  << "\n";
1120  std::cout << " " << diag(i) << "\n";
1121  Index j = (i < k /*|| l==0*/) ? i : perm(l - 1);
1122  std::cout << " "
1123  << "j=" << j << "\n";
1124  }
1125 #endif
1126  Index j = i < k ? i : l > 0 ? perm(l - 1) : i;
1127 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1128  if (!(dk != Literal(0) || diag(i) != Literal(0))) {
1129  std::cout << "k=" << k << ", i=" << i << ", l=" << l << ", perm.size()=" << perm.size() << "\n";
1130  }
1131  eigen_internal_assert(dk != Literal(0) || diag(i) != Literal(0));
1132 #endif
1133  prod *= ((singVals(j) + dk) / ((diag(i) + dk))) * ((mus(j) + (shifts(j) - dk)) / ((diag(i) - dk)));
1134 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1135  eigen_internal_assert(prod >= 0);
1136 #endif
1137 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1138  if (i != k &&
1139  numext::abs(((singVals(j) + dk) * (mus(j) + (shifts(j) - dk))) / ((diag(i) + dk) * (diag(i) - dk)) - 1) >
1140  0.9)
1141  std::cout << " "
1142  << ((singVals(j) + dk) * (mus(j) + (shifts(j) - dk))) / ((diag(i) + dk) * (diag(i) - dk))
1143  << " == (" << (singVals(j) + dk) << " * " << (mus(j) + (shifts(j) - dk)) << ") / ("
1144  << (diag(i) + dk) << " * " << (diag(i) - dk) << ")\n";
1145 #endif
1146  }
1147  }
1148 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1149  std::cout << "zhat(" << k << ") = sqrt( " << prod << ") ; " << (singVals(lastIdx) + dk) << " * "
1150  << mus(lastIdx) + shifts(lastIdx) << " - " << dk << "\n";
1151 #endif
1152  RealScalar tmp = sqrt(prod);
1153 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1154  eigen_internal_assert((numext::isfinite)(tmp));
1155 #endif
1156  zhat(k) = col0(k) > Literal(0) ? RealScalar(tmp) : RealScalar(-tmp);
1157  }
1158  }
1159 }
1160 
1161 // compute singular vectors
1162 template <typename MatrixType, int Options>
1163 void BDCSVD<MatrixType, Options>::computeSingVecs(const ArrayRef& zhat, const ArrayRef& diag, const IndicesRef& perm,
1164  const VectorType& singVals, const ArrayRef& shifts,
1165  const ArrayRef& mus, MatrixXr& U, MatrixXr& V) {
1166  Index n = zhat.size();
1167  Index m = perm.size();
1168 
1169  for (Index k = 0; k < n; ++k) {
1170  if (numext::is_exactly_zero(zhat(k))) {
1171  U.col(k) = VectorType::Unit(n + 1, k);
1172  if (m_compV) V.col(k) = VectorType::Unit(n, k);
1173  } else {
1174  U.col(k).setZero();
1175  for (Index l = 0; l < m; ++l) {
1176  Index i = perm(l);
1177  U(i, k) = zhat(i) / (((diag(i) - shifts(k)) - mus(k))) / ((diag(i) + singVals[k]));
1178  }
1179  U(n, k) = Literal(0);
1180  U.col(k).normalize();
1181 
1182  if (m_compV) {
1183  V.col(k).setZero();
1184  for (Index l = 1; l < m; ++l) {
1185  Index i = perm(l);
1186  V(i, k) = diag(i) * zhat(i) / (((diag(i) - shifts(k)) - mus(k))) / ((diag(i) + singVals[k]));
1187  }
1188  V(0, k) = Literal(-1);
1189  V.col(k).normalize();
1190  }
1191  }
1192  }
1193  U.col(n) = VectorType::Unit(n + 1, n);
1194 }
1195 
1196 // page 12_13
1197 // i >= 1, di almost null and zi non null.
1198 // We use a rotation to zero out zi applied to the left of M, and set di = 0.
1199 template <typename MatrixType, int Options>
1200 void BDCSVD<MatrixType, Options>::deflation43(Index firstCol, Index shift, Index i, Index size) {
1201  using std::abs;
1202  using std::pow;
1203  using std::sqrt;
1204  Index start = firstCol + shift;
1205  RealScalar c = m_computed(start, start);
1206  RealScalar s = m_computed(start + i, start);
1207  RealScalar r = numext::hypot(c, s);
1208  if (numext::is_exactly_zero(r)) {
1209  m_computed(start + i, start + i) = Literal(0);
1210  return;
1211  }
1212  m_computed(start, start) = r;
1213  m_computed(start + i, start) = Literal(0);
1214  m_computed(start + i, start + i) = Literal(0);
1215 
1216  JacobiRotation<RealScalar> J(c / r, -s / r);
1217  if (m_compU)
1218  m_naiveU.middleRows(firstCol, size + 1).applyOnTheRight(firstCol, firstCol + i, J);
1219  else
1220  m_naiveU.applyOnTheRight(firstCol, firstCol + i, J);
1221 } // end deflation 43
1222 
1223 // page 13
1224 // i,j >= 1, i > j, and |di - dj| < epsilon * norm2(M)
1225 // We apply two rotations to have zi = 0, and dj = di.
1226 template <typename MatrixType, int Options>
1227 void BDCSVD<MatrixType, Options>::deflation44(Index firstColu, Index firstColm, Index firstRowW, Index firstColW,
1228  Index i, Index j, Index size) {
1229  using std::abs;
1230  using std::conj;
1231  using std::pow;
1232  using std::sqrt;
1233 
1234  RealScalar s = m_computed(firstColm + i, firstColm);
1235  RealScalar c = m_computed(firstColm + j, firstColm);
1236  RealScalar r = numext::hypot(c, s);
1237 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1238  std::cout << "deflation 4.4: " << i << "," << j << " -> " << c << " " << s << " " << r << " ; "
1239  << m_computed(firstColm + i - 1, firstColm) << " " << m_computed(firstColm + i, firstColm) << " "
1240  << m_computed(firstColm + i + 1, firstColm) << " " << m_computed(firstColm + i + 2, firstColm) << "\n";
1241  std::cout << m_computed(firstColm + i - 1, firstColm + i - 1) << " " << m_computed(firstColm + i, firstColm + i)
1242  << " " << m_computed(firstColm + i + 1, firstColm + i + 1) << " "
1243  << m_computed(firstColm + i + 2, firstColm + i + 2) << "\n";
1244 #endif
1245  if (numext::is_exactly_zero(r)) {
1246  m_computed(firstColm + j, firstColm + j) = m_computed(firstColm + i, firstColm + i);
1247  return;
1248  }
1249  c /= r;
1250  s /= r;
1251  m_computed(firstColm + j, firstColm) = r;
1252  m_computed(firstColm + j, firstColm + j) = m_computed(firstColm + i, firstColm + i);
1253  m_computed(firstColm + i, firstColm) = Literal(0);
1254 
1255  JacobiRotation<RealScalar> J(c, -s);
1256  if (m_compU)
1257  m_naiveU.middleRows(firstColu, size + 1).applyOnTheRight(firstColu + j, firstColu + i, J);
1258  else
1259  m_naiveU.applyOnTheRight(firstColu + j, firstColu + i, J);
1260  if (m_compV) m_naiveV.middleRows(firstRowW, size).applyOnTheRight(firstColW + j, firstColW + i, J);
1261 } // end deflation 44
1262 
1263 // acts on block from (firstCol+shift, firstCol+shift) to (lastCol+shift, lastCol+shift) [inclusive]
1264 template <typename MatrixType, int Options>
1265 void BDCSVD<MatrixType, Options>::deflation(Index firstCol, Index lastCol, Index k, Index firstRowW, Index firstColW,
1266  Index shift) {
1267  using std::abs;
1268  using std::sqrt;
1269  const Index length = lastCol + 1 - firstCol;
1270 
1271  Block<MatrixXr, Dynamic, 1> col0(m_computed, firstCol + shift, firstCol + shift, length, 1);
1272  Diagonal<MatrixXr> fulldiag(m_computed);
1273  VectorBlock<Diagonal<MatrixXr>, Dynamic> diag(fulldiag, firstCol + shift, length);
1274 
1275  const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
1276  RealScalar maxDiag = diag.tail((std::max)(Index(1), length - 1)).cwiseAbs().maxCoeff();
1277  RealScalar epsilon_strict = numext::maxi<RealScalar>(considerZero, NumTraits<RealScalar>::epsilon() * maxDiag);
1278  RealScalar epsilon_coarse =
1279  Literal(8) * NumTraits<RealScalar>::epsilon() * numext::maxi<RealScalar>(col0.cwiseAbs().maxCoeff(), maxDiag);
1280 
1281 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1282  eigen_internal_assert(m_naiveU.allFinite());
1283  eigen_internal_assert(m_naiveV.allFinite());
1284  eigen_internal_assert(m_computed.allFinite());
1285 #endif
1286 
1287 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1288  std::cout << "\ndeflate:" << diag.head(k + 1).transpose() << " | "
1289  << diag.segment(k + 1, length - k - 1).transpose() << "\n";
1290 #endif
1291 
1292  // condition 4.1
1293  if (diag(0) < epsilon_coarse) {
1294 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1295  std::cout << "deflation 4.1, because " << diag(0) << " < " << epsilon_coarse << "\n";
1296 #endif
1297  diag(0) = epsilon_coarse;
1298  }
1299 
1300  // condition 4.2
1301  for (Index i = 1; i < length; ++i)
1302  if (abs(col0(i)) < epsilon_strict) {
1303 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1304  std::cout << "deflation 4.2, set z(" << i << ") to zero because " << abs(col0(i)) << " < " << epsilon_strict
1305  << " (diag(" << i << ")=" << diag(i) << ")\n";
1306 #endif
1307  col0(i) = Literal(0);
1308  }
1309 
1310  // condition 4.3
1311  for (Index i = 1; i < length; i++)
1312  if (diag(i) < epsilon_coarse) {
1313 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1314  std::cout << "deflation 4.3, cancel z(" << i << ")=" << col0(i) << " because diag(" << i << ")=" << diag(i)
1315  << " < " << epsilon_coarse << "\n";
1316 #endif
1317  deflation43(firstCol, shift, i, length);
1318  }
1319 
1320 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1321  eigen_internal_assert(m_naiveU.allFinite());
1322  eigen_internal_assert(m_naiveV.allFinite());
1323  eigen_internal_assert(m_computed.allFinite());
1324 #endif
1325 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1326  std::cout << "to be sorted: " << diag.transpose() << "\n\n";
1327  std::cout << " : " << col0.transpose() << "\n\n";
1328 #endif
1329  {
1330  // Check for total deflation:
1331  // If we have a total deflation, then we have to consider col0(0)==diag(0) as a singular value during sorting.
1332  const bool total_deflation = (col0.tail(length - 1).array().abs() < considerZero).all();
1333 
1334  // Sort the diagonal entries, since diag(1:k-1) and diag(k:length) are already sorted, let's do a sorted merge.
1335  // First, compute the respective permutation.
1336  Index* permutation = m_workspaceI.data();
1337  {
1338  permutation[0] = 0;
1339  Index p = 1;
1340 
1341  // Move deflated diagonal entries at the end.
1342  for (Index i = 1; i < length; ++i)
1343  if (diag(i) < considerZero) permutation[p++] = i;
1344 
1345  Index i = 1, j = k + 1;
1346  for (; p < length; ++p) {
1347  if (i > k)
1348  permutation[p] = j++;
1349  else if (j >= length)
1350  permutation[p] = i++;
1351  else if (diag(i) < diag(j))
1352  permutation[p] = j++;
1353  else
1354  permutation[p] = i++;
1355  }
1356  }
1357 
1358  // If we have a total deflation, then we have to insert diag(0) at the right place
1359  if (total_deflation) {
1360  for (Index i = 1; i < length; ++i) {
1361  Index pi = permutation[i];
1362  if (diag(pi) < considerZero || diag(0) < diag(pi))
1363  permutation[i - 1] = permutation[i];
1364  else {
1365  permutation[i - 1] = 0;
1366  break;
1367  }
1368  }
1369  }
1370 
1371  // Current index of each col, and current column of each index
1372  Index* realInd = m_workspaceI.data() + length;
1373  Index* realCol = m_workspaceI.data() + 2 * length;
1374 
1375  for (int pos = 0; pos < length; pos++) {
1376  realCol[pos] = pos;
1377  realInd[pos] = pos;
1378  }
1379 
1380  for (Index i = total_deflation ? 0 : 1; i < length; i++) {
1381  const Index pi = permutation[length - (total_deflation ? i + 1 : i)];
1382  const Index J = realCol[pi];
1383 
1384  using std::swap;
1385  // swap diagonal and first column entries:
1386  swap(diag(i), diag(J));
1387  if (i != 0 && J != 0) swap(col0(i), col0(J));
1388 
1389  // change columns
1390  if (m_compU)
1391  m_naiveU.col(firstCol + i)
1392  .segment(firstCol, length + 1)
1393  .swap(m_naiveU.col(firstCol + J).segment(firstCol, length + 1));
1394  else
1395  m_naiveU.col(firstCol + i).segment(0, 2).swap(m_naiveU.col(firstCol + J).segment(0, 2));
1396  if (m_compV)
1397  m_naiveV.col(firstColW + i)
1398  .segment(firstRowW, length)
1399  .swap(m_naiveV.col(firstColW + J).segment(firstRowW, length));
1400 
1401  // update real pos
1402  const Index realI = realInd[i];
1403  realCol[realI] = J;
1404  realCol[pi] = i;
1405  realInd[J] = realI;
1406  realInd[i] = pi;
1407  }
1408  }
1409 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1410  std::cout << "sorted: " << diag.transpose().format(bdcsvdfmt) << "\n";
1411  std::cout << " : " << col0.transpose() << "\n\n";
1412 #endif
1413 
1414  // condition 4.4
1415  {
1416  Index i = length - 1;
1417  // Find last non-deflated entry.
1418  while (i > 0 && (diag(i) < considerZero || abs(col0(i)) < considerZero)) --i;
1419 
1420  for (; i > 1; --i)
1421  if ((diag(i) - diag(i - 1)) < epsilon_strict) {
1422 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1423  std::cout << "deflation 4.4 with i = " << i << " because " << diag(i) << " - " << diag(i - 1)
1424  << " == " << (diag(i) - diag(i - 1)) << " < " << epsilon_strict << "\n";
1425 #endif
1426  eigen_internal_assert(abs(diag(i) - diag(i - 1)) < epsilon_coarse &&
1427  " diagonal entries are not properly sorted");
1428  deflation44(firstCol, firstCol + shift, firstRowW, firstColW, i, i - 1, length);
1429  }
1430  }
1431 
1432 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1433  for (Index j = 2; j < length; ++j) eigen_internal_assert(diag(j - 1) <= diag(j) || abs(diag(j)) < considerZero);
1434 #endif
1435 
1436 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1437  eigen_internal_assert(m_naiveU.allFinite());
1438  eigen_internal_assert(m_naiveV.allFinite());
1439  eigen_internal_assert(m_computed.allFinite());
1440 #endif
1441 } // end deflation
1442 
1449 template <typename Derived>
1450 template <int Options>
1452  return BDCSVD<PlainObject, Options>(*this);
1453 }
1454 
1461 template <typename Derived>
1462 template <int Options>
1464  unsigned int computationOptions) const {
1465  return BDCSVD<PlainObject, Options>(*this, computationOptions);
1466 }
1467 
1468 } // end namespace Eigen
1469 
1470 #endif
static constexpr Eigen::internal::all_t all
Definition: IndexedViewHelper.h:86
constexpr Index size() const noexcept
Definition: EigenBase.h:64
Matrix< Type, Dynamic, Dynamic > MatrixX
[c++11] Dynamic×Dynamic matrix of type Type.
Definition: Matrix.h:514
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
BDCSVD & compute(const MatrixBase< Derived > &matrix)
Method performing the decomposition of given matrix. Computes Thin/Full unitaries U/V if specified us...
Definition: BDCSVD.h:201
bool computeV() const
Definition: SVDBase.h:275
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
std::enable_if_t< std::is_base_of< DenseBase< std::decay_t< DerivedA > >, std::decay_t< DerivedA > >::value &&std::is_base_of< DenseBase< std::decay_t< DerivedB > >, std::decay_t< DerivedB > >::value, void > swap(DerivedA &&a, DerivedB &&b)
Definition: DenseBase.h:667
Definition: Constants.h:242
Definition: Constants.h:447
class Bidiagonal Divide and Conquer SVD
Definition: ForwardDeclarations.h:439
Definition: Constants.h:440
EIGEN_DEPRECATED_WITH_REASON("Initialization is no longer needed.") inline void initParallel()
Definition: Parallelizer.h:50
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
BDCSVD(const MatrixBase< Derived > &matrix)
Constructor performing the decomposition of given matrix, using the custom options specified with the...
Definition: BDCSVD.h:170
bool computeU() const
Definition: SVDBase.h:273
Definition: Constants.h:429
const int Dynamic
Definition: Constants.h:25
BDCSVD(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: BDCSVD.h:138
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
Definition: Constants.h:444
BDCSVD()
Default Constructor.
Definition: BDCSVD.h:130