• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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-2016 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 namespace Eigen {
26 
27 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
28 IOFormat bdcsvdfmt(8, 0, ", ", "\n", "  [", "]");
29 #endif
30 
31 template<typename _MatrixType> class BDCSVD;
32 
33 namespace internal {
34 
35 template<typename _MatrixType>
36 struct traits<BDCSVD<_MatrixType> >
37 {
38   typedef _MatrixType MatrixType;
39 };
40 
41 } // end namespace internal
42 
43 
44 /** \ingroup SVD_Module
45  *
46  *
47  * \class BDCSVD
48  *
49  * \brief class Bidiagonal Divide and Conquer SVD
50  *
51  * \tparam _MatrixType the type of the matrix of which we are computing the SVD decomposition
52  *
53  * This class first reduces the input matrix to bi-diagonal form using class UpperBidiagonalization,
54  * and then performs a divide-and-conquer diagonalization. Small blocks are diagonalized using class JacobiSVD.
55  * You can control the switching size with the setSwitchSize() method, default is 16.
56  * For small matrice (<16), it is thus preferable to directly use JacobiSVD. For larger ones, BDCSVD is highly
57  * recommended and can several order of magnitude faster.
58  *
59  * \warning this algorithm is unlikely to provide accurate result when compiled with unsafe math optimizations.
60  * For instance, this concerns Intel's compiler (ICC), which perfroms such optimization by default unless
61  * you compile with the \c -fp-model \c precise option. Likewise, the \c -ffast-math option of GCC or clang will
62  * significantly degrade the accuracy.
63  *
64  * \sa class JacobiSVD
65  */
66 template<typename _MatrixType>
67 class BDCSVD : public SVDBase<BDCSVD<_MatrixType> >
68 {
69   typedef SVDBase<BDCSVD> Base;
70 
71 public:
72   using Base::rows;
73   using Base::cols;
74   using Base::computeU;
75   using Base::computeV;
76 
77   typedef _MatrixType MatrixType;
78   typedef typename MatrixType::Scalar Scalar;
79   typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
80   enum {
81     RowsAtCompileTime = MatrixType::RowsAtCompileTime,
82     ColsAtCompileTime = MatrixType::ColsAtCompileTime,
83     DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime, ColsAtCompileTime),
84     MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
85     MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
86     MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime, MaxColsAtCompileTime),
87     MatrixOptions = MatrixType::Options
88   };
89 
90   typedef typename Base::MatrixUType MatrixUType;
91   typedef typename Base::MatrixVType MatrixVType;
92   typedef typename Base::SingularValuesType SingularValuesType;
93 
94   typedef Matrix<Scalar, Dynamic, Dynamic, ColMajor> MatrixX;
95   typedef Matrix<RealScalar, Dynamic, Dynamic, ColMajor> MatrixXr;
96   typedef Matrix<RealScalar, Dynamic, 1> VectorType;
97   typedef Array<RealScalar, Dynamic, 1> ArrayXr;
98   typedef Array<Index,1,Dynamic> ArrayXi;
99   typedef Ref<ArrayXr> ArrayRef;
100   typedef Ref<ArrayXi> IndicesRef;
101 
102   /** \brief Default Constructor.
103    *
104    * The default constructor is useful in cases in which the user intends to
105    * perform decompositions via BDCSVD::compute(const MatrixType&).
106    */
107   BDCSVD() : m_algoswap(16), m_numIters(0)
108   {}
109 
110 
111   /** \brief Default Constructor with memory preallocation
112    *
113    * Like the default constructor but with preallocation of the internal data
114    * according to the specified problem size.
115    * \sa BDCSVD()
116    */
117   BDCSVD(Index rows, Index cols, unsigned int computationOptions = 0)
118     : m_algoswap(16), m_numIters(0)
119   {
120     allocate(rows, cols, computationOptions);
121   }
122 
123   /** \brief Constructor performing the decomposition of given matrix.
124    *
125    * \param matrix the matrix to decompose
126    * \param computationOptions optional parameter allowing to specify if you want full or thin U or V unitaries to be computed.
127    *                           By default, none is computed. This is a bit - field, the possible bits are #ComputeFullU, #ComputeThinU,
128    *                           #ComputeFullV, #ComputeThinV.
129    *
130    * Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not
131    * available with the (non - default) FullPivHouseholderQR preconditioner.
132    */
133   BDCSVD(const MatrixType& matrix, unsigned int computationOptions = 0)
134     : m_algoswap(16), m_numIters(0)
135   {
136     compute(matrix, computationOptions);
137   }
138 
139   ~BDCSVD()
140   {
141   }
142 
143   /** \brief Method performing the decomposition of given matrix using custom options.
144    *
145    * \param matrix the matrix to decompose
146    * \param computationOptions optional parameter allowing to specify if you want full or thin U or V unitaries to be computed.
147    *                           By default, none is computed. This is a bit - field, the possible bits are #ComputeFullU, #ComputeThinU,
148    *                           #ComputeFullV, #ComputeThinV.
149    *
150    * Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not
151    * available with the (non - default) FullPivHouseholderQR preconditioner.
152    */
153   BDCSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
154 
155   /** \brief Method performing the decomposition of given matrix using current options.
156    *
157    * \param matrix the matrix to decompose
158    *
159    * This method uses the current \a computationOptions, as already passed to the constructor or to compute(const MatrixType&, unsigned int).
160    */
161   BDCSVD& compute(const MatrixType& matrix)
162   {
163     return compute(matrix, this->m_computationOptions);
164   }
165 
166   void setSwitchSize(int s)
167   {
168     eigen_assert(s>3 && "BDCSVD the size of the algo switch has to be greater than 3");
169     m_algoswap = s;
170   }
171 
172 private:
173   void allocate(Index rows, Index cols, unsigned int computationOptions);
174   void divide(Index firstCol, Index lastCol, Index firstRowW, Index firstColW, Index shift);
175   void computeSVDofM(Index firstCol, Index n, MatrixXr& U, VectorType& singVals, MatrixXr& V);
176   void computeSingVals(const ArrayRef& col0, const ArrayRef& diag, const IndicesRef& perm, VectorType& singVals, ArrayRef shifts, ArrayRef mus);
177   void perturbCol0(const ArrayRef& col0, const ArrayRef& diag, const IndicesRef& perm, const VectorType& singVals, const ArrayRef& shifts, const ArrayRef& mus, ArrayRef zhat);
178   void computeSingVecs(const ArrayRef& zhat, const ArrayRef& diag, const IndicesRef& perm, const VectorType& singVals, const ArrayRef& shifts, const ArrayRef& mus, MatrixXr& U, MatrixXr& V);
179   void deflation43(Index firstCol, Index shift, Index i, Index size);
180   void deflation44(Index firstColu , Index firstColm, Index firstRowW, Index firstColW, Index i, Index j, Index size);
181   void deflation(Index firstCol, Index lastCol, Index k, Index firstRowW, Index firstColW, Index shift);
182   template<typename HouseholderU, typename HouseholderV, typename NaiveU, typename NaiveV>
183   void copyUV(const HouseholderU &householderU, const HouseholderV &householderV, const NaiveU &naiveU, const NaiveV &naivev);
184   void structured_update(Block<MatrixXr,Dynamic,Dynamic> A, const MatrixXr &B, Index n1);
185   static RealScalar secularEq(RealScalar x, const ArrayRef& col0, const ArrayRef& diag, const IndicesRef &perm, const ArrayRef& diagShifted, RealScalar shift);
186 
187 protected:
188   MatrixXr m_naiveU, m_naiveV;
189   MatrixXr m_computed;
190   Index m_nRec;
191   ArrayXr m_workspace;
192   ArrayXi m_workspaceI;
193   int m_algoswap;
194   bool m_isTranspose, m_compU, m_compV;
195 
196   using Base::m_singularValues;
197   using Base::m_diagSize;
198   using Base::m_computeFullU;
199   using Base::m_computeFullV;
200   using Base::m_computeThinU;
201   using Base::m_computeThinV;
202   using Base::m_matrixU;
203   using Base::m_matrixV;
204   using Base::m_isInitialized;
205   using Base::m_nonzeroSingularValues;
206 
207 public:
208   int m_numIters;
209 }; //end class BDCSVD
210 
211 
212 // Method to allocate and initialize matrix and attributes
213 template<typename MatrixType>
214 void BDCSVD<MatrixType>::allocate(Index rows, Index cols, unsigned int computationOptions)
215 {
216   m_isTranspose = (cols > rows);
217 
218   if (Base::allocate(rows, cols, computationOptions))
219     return;
220 
221   m_computed = MatrixXr::Zero(m_diagSize + 1, m_diagSize );
222   m_compU = computeV();
223   m_compV = computeU();
224   if (m_isTranspose)
225     std::swap(m_compU, m_compV);
226 
227   if (m_compU) m_naiveU = MatrixXr::Zero(m_diagSize + 1, m_diagSize + 1 );
228   else         m_naiveU = MatrixXr::Zero(2, m_diagSize + 1 );
229 
230   if (m_compV) m_naiveV = MatrixXr::Zero(m_diagSize, m_diagSize);
231 
232   m_workspace.resize((m_diagSize+1)*(m_diagSize+1)*3);
233   m_workspaceI.resize(3*m_diagSize);
234 }// end allocate
235 
236 template<typename MatrixType>
237 BDCSVD<MatrixType>& BDCSVD<MatrixType>::compute(const MatrixType& matrix, unsigned int computationOptions)
238 {
239 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
240   std::cout << "\n\n\n======================================================================================================================\n\n\n";
241 #endif
242   allocate(matrix.rows(), matrix.cols(), computationOptions);
243   using std::abs;
244 
245   const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
246 
247   //**** step -1 - If the problem is too small, directly falls back to JacobiSVD and return
248   if(matrix.cols() < m_algoswap)
249   {
250     // FIXME this line involves temporaries
251     JacobiSVD<MatrixType> jsvd(matrix,computationOptions);
252     if(computeU()) m_matrixU = jsvd.matrixU();
253     if(computeV()) m_matrixV = jsvd.matrixV();
254     m_singularValues = jsvd.singularValues();
255     m_nonzeroSingularValues = jsvd.nonzeroSingularValues();
256     m_isInitialized = true;
257     return *this;
258   }
259 
260   //**** step 0 - Copy the input matrix and apply scaling to reduce over/under-flows
261   RealScalar scale = matrix.cwiseAbs().maxCoeff();
262   if(scale==RealScalar(0)) scale = RealScalar(1);
263   MatrixX copy;
264   if (m_isTranspose) copy = matrix.adjoint()/scale;
265   else               copy = matrix/scale;
266 
267   //**** step 1 - Bidiagonalization
268   // FIXME this line involves temporaries
269   internal::UpperBidiagonalization<MatrixX> bid(copy);
270 
271   //**** step 2 - Divide & Conquer
272   m_naiveU.setZero();
273   m_naiveV.setZero();
274   // FIXME this line involves a temporary matrix
275   m_computed.topRows(m_diagSize) = bid.bidiagonal().toDenseMatrix().transpose();
276   m_computed.template bottomRows<1>().setZero();
277   divide(0, m_diagSize - 1, 0, 0, 0);
278 
279   //**** step 3 - Copy singular values and vectors
280   for (int i=0; i<m_diagSize; i++)
281   {
282     RealScalar a = abs(m_computed.coeff(i, i));
283     m_singularValues.coeffRef(i) = a * scale;
284     if (a<considerZero)
285     {
286       m_nonzeroSingularValues = i;
287       m_singularValues.tail(m_diagSize - i - 1).setZero();
288       break;
289     }
290     else if (i == m_diagSize - 1)
291     {
292       m_nonzeroSingularValues = i + 1;
293       break;
294     }
295   }
296 
297 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
298 //   std::cout << "m_naiveU\n" << m_naiveU << "\n\n";
299 //   std::cout << "m_naiveV\n" << m_naiveV << "\n\n";
300 #endif
301   if(m_isTranspose) copyUV(bid.householderV(), bid.householderU(), m_naiveV, m_naiveU);
302   else              copyUV(bid.householderU(), bid.householderV(), m_naiveU, m_naiveV);
303 
304   m_isInitialized = true;
305   return *this;
306 }// end compute
307 
308 
309 template<typename MatrixType>
310 template<typename HouseholderU, typename HouseholderV, typename NaiveU, typename NaiveV>
311 void BDCSVD<MatrixType>::copyUV(const HouseholderU &householderU, const HouseholderV &householderV, const NaiveU &naiveU, const NaiveV &naiveV)
312 {
313   // Note exchange of U and V: m_matrixU is set from m_naiveV and vice versa
314   if (computeU())
315   {
316     Index Ucols = m_computeThinU ? m_diagSize : householderU.cols();
317     m_matrixU = MatrixX::Identity(householderU.cols(), Ucols);
318     m_matrixU.topLeftCorner(m_diagSize, m_diagSize) = naiveV.template cast<Scalar>().topLeftCorner(m_diagSize, m_diagSize);
319     householderU.applyThisOnTheLeft(m_matrixU); // FIXME this line involves a temporary buffer
320   }
321   if (computeV())
322   {
323     Index Vcols = m_computeThinV ? m_diagSize : householderV.cols();
324     m_matrixV = MatrixX::Identity(householderV.cols(), Vcols);
325     m_matrixV.topLeftCorner(m_diagSize, m_diagSize) = naiveU.template cast<Scalar>().topLeftCorner(m_diagSize, m_diagSize);
326     householderV.applyThisOnTheLeft(m_matrixV); // FIXME this line involves a temporary buffer
327   }
328 }
329 
330 /** \internal
331   * Performs A = A * B exploiting the special structure of the matrix A. Splitting A as:
332   *  A = [A1]
333   *      [A2]
334   * such that A1.rows()==n1, then we assume that at least half of the columns of A1 and A2 are zeros.
335   * We can thus pack them prior to the the matrix product. However, this is only worth the effort if the matrix is large
336   * enough.
337   */
338 template<typename MatrixType>
339 void BDCSVD<MatrixType>::structured_update(Block<MatrixXr,Dynamic,Dynamic> A, const MatrixXr &B, Index n1)
340 {
341   Index n = A.rows();
342   if(n>100)
343   {
344     // If the matrices are large enough, let's exploit the sparse structure of A by
345     // splitting it in half (wrt n1), and packing the non-zero columns.
346     Index n2 = n - n1;
347     Map<MatrixXr> A1(m_workspace.data()      , n1, n);
348     Map<MatrixXr> A2(m_workspace.data()+ n1*n, n2, n);
349     Map<MatrixXr> B1(m_workspace.data()+  n*n, n,  n);
350     Map<MatrixXr> B2(m_workspace.data()+2*n*n, n,  n);
351     Index k1=0, k2=0;
352     for(Index j=0; j<n; ++j)
353     {
354       if( (A.col(j).head(n1).array()!=0).any() )
355       {
356         A1.col(k1) = A.col(j).head(n1);
357         B1.row(k1) = B.row(j);
358         ++k1;
359       }
360       if( (A.col(j).tail(n2).array()!=0).any() )
361       {
362         A2.col(k2) = A.col(j).tail(n2);
363         B2.row(k2) = B.row(j);
364         ++k2;
365       }
366     }
367 
368     A.topRows(n1).noalias()    = A1.leftCols(k1) * B1.topRows(k1);
369     A.bottomRows(n2).noalias() = A2.leftCols(k2) * B2.topRows(k2);
370   }
371   else
372   {
373     Map<MatrixXr,Aligned> tmp(m_workspace.data(),n,n);
374     tmp.noalias() = A*B;
375     A = tmp;
376   }
377 }
378 
379 // The divide algorithm is done "in place", we are always working on subsets of the same matrix. The divide methods takes as argument the
380 // place of the submatrix we are currently working on.
381 
382 //@param firstCol : The Index of the first column of the submatrix of m_computed and for m_naiveU;
383 //@param lastCol : The Index of the last column of the submatrix of m_computed and for m_naiveU;
384 // lastCol + 1 - firstCol is the size of the submatrix.
385 //@param firstRowW : The Index of the first row of the matrix W that we are to change. (see the reference paper section 1 for more information on W)
386 //@param firstRowW : Same as firstRowW with the column.
387 //@param shift : Each time one takes the left submatrix, one must add 1 to the shift. Why? Because! We actually want the last column of the U submatrix
388 // to become the first column (*coeff) and to shift all the other columns to the right. There are more details on the reference paper.
389 template<typename MatrixType>
390 void BDCSVD<MatrixType>::divide (Index firstCol, Index lastCol, Index firstRowW, Index firstColW, Index shift)
391 {
392   // requires rows = cols + 1;
393   using std::pow;
394   using std::sqrt;
395   using std::abs;
396   const Index n = lastCol - firstCol + 1;
397   const Index k = n/2;
398   const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
399   RealScalar alphaK;
400   RealScalar betaK;
401   RealScalar r0;
402   RealScalar lambda, phi, c0, s0;
403   VectorType l, f;
404   // We use the other algorithm which is more efficient for small
405   // matrices.
406   if (n < m_algoswap)
407   {
408     // FIXME this line involves temporaries
409     JacobiSVD<MatrixXr> b(m_computed.block(firstCol, firstCol, n + 1, n), ComputeFullU | (m_compV ? ComputeFullV : 0));
410     if (m_compU)
411       m_naiveU.block(firstCol, firstCol, n + 1, n + 1).real() = b.matrixU();
412     else
413     {
414       m_naiveU.row(0).segment(firstCol, n + 1).real() = b.matrixU().row(0);
415       m_naiveU.row(1).segment(firstCol, n + 1).real() = b.matrixU().row(n);
416     }
417     if (m_compV) m_naiveV.block(firstRowW, firstColW, n, n).real() = b.matrixV();
418     m_computed.block(firstCol + shift, firstCol + shift, n + 1, n).setZero();
419     m_computed.diagonal().segment(firstCol + shift, n) = b.singularValues().head(n);
420     return;
421   }
422   // We use the divide and conquer algorithm
423   alphaK =  m_computed(firstCol + k, firstCol + k);
424   betaK = m_computed(firstCol + k + 1, firstCol + k);
425   // The divide must be done in that order in order to have good results. Divide change the data inside the submatrices
426   // and the divide of the right submatrice reads one column of the left submatrice. That's why we need to treat the
427   // right submatrix before the left one.
428   divide(k + 1 + firstCol, lastCol, k + 1 + firstRowW, k + 1 + firstColW, shift);
429   divide(firstCol, k - 1 + firstCol, firstRowW, firstColW + 1, shift + 1);
430 
431   if (m_compU)
432   {
433     lambda = m_naiveU(firstCol + k, firstCol + k);
434     phi = m_naiveU(firstCol + k + 1, lastCol + 1);
435   }
436   else
437   {
438     lambda = m_naiveU(1, firstCol + k);
439     phi = m_naiveU(0, lastCol + 1);
440   }
441   r0 = sqrt((abs(alphaK * lambda) * abs(alphaK * lambda)) + abs(betaK * phi) * abs(betaK * phi));
442   if (m_compU)
443   {
444     l = m_naiveU.row(firstCol + k).segment(firstCol, k);
445     f = m_naiveU.row(firstCol + k + 1).segment(firstCol + k + 1, n - k - 1);
446   }
447   else
448   {
449     l = m_naiveU.row(1).segment(firstCol, k);
450     f = m_naiveU.row(0).segment(firstCol + k + 1, n - k - 1);
451   }
452   if (m_compV) m_naiveV(firstRowW+k, firstColW) = 1;
453   if (r0<considerZero)
454   {
455     c0 = 1;
456     s0 = 0;
457   }
458   else
459   {
460     c0 = alphaK * lambda / r0;
461     s0 = betaK * phi / r0;
462   }
463 
464 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
465   assert(m_naiveU.allFinite());
466   assert(m_naiveV.allFinite());
467   assert(m_computed.allFinite());
468 #endif
469 
470   if (m_compU)
471   {
472     MatrixXr q1 (m_naiveU.col(firstCol + k).segment(firstCol, k + 1));
473     // we shiftW Q1 to the right
474     for (Index i = firstCol + k - 1; i >= firstCol; i--)
475       m_naiveU.col(i + 1).segment(firstCol, k + 1) = m_naiveU.col(i).segment(firstCol, k + 1);
476     // we shift q1 at the left with a factor c0
477     m_naiveU.col(firstCol).segment( firstCol, k + 1) = (q1 * c0);
478     // last column = q1 * - s0
479     m_naiveU.col(lastCol + 1).segment(firstCol, k + 1) = (q1 * ( - s0));
480     // first column = q2 * s0
481     m_naiveU.col(firstCol).segment(firstCol + k + 1, n - k) = m_naiveU.col(lastCol + 1).segment(firstCol + k + 1, n - k) * s0;
482     // q2 *= c0
483     m_naiveU.col(lastCol + 1).segment(firstCol + k + 1, n - k) *= c0;
484   }
485   else
486   {
487     RealScalar q1 = m_naiveU(0, firstCol + k);
488     // we shift Q1 to the right
489     for (Index i = firstCol + k - 1; i >= firstCol; i--)
490       m_naiveU(0, i + 1) = m_naiveU(0, i);
491     // we shift q1 at the left with a factor c0
492     m_naiveU(0, firstCol) = (q1 * c0);
493     // last column = q1 * - s0
494     m_naiveU(0, lastCol + 1) = (q1 * ( - s0));
495     // first column = q2 * s0
496     m_naiveU(1, firstCol) = m_naiveU(1, lastCol + 1) *s0;
497     // q2 *= c0
498     m_naiveU(1, lastCol + 1) *= c0;
499     m_naiveU.row(1).segment(firstCol + 1, k).setZero();
500     m_naiveU.row(0).segment(firstCol + k + 1, n - k - 1).setZero();
501   }
502 
503 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
504   assert(m_naiveU.allFinite());
505   assert(m_naiveV.allFinite());
506   assert(m_computed.allFinite());
507 #endif
508 
509   m_computed(firstCol + shift, firstCol + shift) = r0;
510   m_computed.col(firstCol + shift).segment(firstCol + shift + 1, k) = alphaK * l.transpose().real();
511   m_computed.col(firstCol + shift).segment(firstCol + shift + k + 1, n - k - 1) = betaK * f.transpose().real();
512 
513 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
514   ArrayXr tmp1 = (m_computed.block(firstCol+shift, firstCol+shift, n, n)).jacobiSvd().singularValues();
515 #endif
516   // Second part: try to deflate singular values in combined matrix
517   deflation(firstCol, lastCol, k, firstRowW, firstColW, shift);
518 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
519   ArrayXr tmp2 = (m_computed.block(firstCol+shift, firstCol+shift, n, n)).jacobiSvd().singularValues();
520   std::cout << "\n\nj1 = " << tmp1.transpose().format(bdcsvdfmt) << "\n";
521   std::cout << "j2 = " << tmp2.transpose().format(bdcsvdfmt) << "\n\n";
522   std::cout << "err:      " << ((tmp1-tmp2).abs()>1e-12*tmp2.abs()).transpose() << "\n";
523   static int count = 0;
524   std::cout << "# " << ++count << "\n\n";
525   assert((tmp1-tmp2).matrix().norm() < 1e-14*tmp2.matrix().norm());
526 //   assert(count<681);
527 //   assert(((tmp1-tmp2).abs()<1e-13*tmp2.abs()).all());
528 #endif
529 
530   // Third part: compute SVD of combined matrix
531   MatrixXr UofSVD, VofSVD;
532   VectorType singVals;
533   computeSVDofM(firstCol + shift, n, UofSVD, singVals, VofSVD);
534 
535 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
536   assert(UofSVD.allFinite());
537   assert(VofSVD.allFinite());
538 #endif
539 
540   if (m_compU)
541     structured_update(m_naiveU.block(firstCol, firstCol, n + 1, n + 1), UofSVD, (n+2)/2);
542   else
543   {
544     Map<Matrix<RealScalar,2,Dynamic>,Aligned> tmp(m_workspace.data(),2,n+1);
545     tmp.noalias() = m_naiveU.middleCols(firstCol, n+1) * UofSVD;
546     m_naiveU.middleCols(firstCol, n + 1) = tmp;
547   }
548 
549   if (m_compV)  structured_update(m_naiveV.block(firstRowW, firstColW, n, n), VofSVD, (n+1)/2);
550 
551 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
552   assert(m_naiveU.allFinite());
553   assert(m_naiveV.allFinite());
554   assert(m_computed.allFinite());
555 #endif
556 
557   m_computed.block(firstCol + shift, firstCol + shift, n, n).setZero();
558   m_computed.block(firstCol + shift, firstCol + shift, n, n).diagonal() = singVals;
559 }// end divide
560 
561 // Compute SVD of m_computed.block(firstCol, firstCol, n + 1, n); this block only has non-zeros in
562 // the first column and on the diagonal and has undergone deflation, so diagonal is in increasing
563 // order except for possibly the (0,0) entry. The computed SVD is stored U, singVals and V, except
564 // that if m_compV is false, then V is not computed. Singular values are sorted in decreasing order.
565 //
566 // TODO Opportunities for optimization: better root finding algo, better stopping criterion, better
567 // handling of round-off errors, be consistent in ordering
568 // For instance, to solve the secular equation using FMM, see http://www.stat.uchicago.edu/~lekheng/courses/302/classics/greengard-rokhlin.pdf
569 template <typename MatrixType>
570 void BDCSVD<MatrixType>::computeSVDofM(Index firstCol, Index n, MatrixXr& U, VectorType& singVals, MatrixXr& V)
571 {
572   const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
573   using std::abs;
574   ArrayRef col0 = m_computed.col(firstCol).segment(firstCol, n);
575   m_workspace.head(n) =  m_computed.block(firstCol, firstCol, n, n).diagonal();
576   ArrayRef diag = m_workspace.head(n);
577   diag(0) = 0;
578 
579   // Allocate space for singular values and vectors
580   singVals.resize(n);
581   U.resize(n+1, n+1);
582   if (m_compV) V.resize(n, n);
583 
584 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
585   if (col0.hasNaN() || diag.hasNaN())
586     std::cout << "\n\nHAS NAN\n\n";
587 #endif
588 
589   // Many singular values might have been deflated, the zero ones have been moved to the end,
590   // but others are interleaved and we must ignore them at this stage.
591   // To this end, let's compute a permutation skipping them:
592   Index actual_n = n;
593   while(actual_n>1 && diag(actual_n-1)==0) --actual_n;
594   Index m = 0; // size of the deflated problem
595   for(Index k=0;k<actual_n;++k)
596     if(abs(col0(k))>considerZero)
597       m_workspaceI(m++) = k;
598   Map<ArrayXi> perm(m_workspaceI.data(),m);
599 
600   Map<ArrayXr> shifts(m_workspace.data()+1*n, n);
601   Map<ArrayXr> mus(m_workspace.data()+2*n, n);
602   Map<ArrayXr> zhat(m_workspace.data()+3*n, n);
603 
604 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
605   std::cout << "computeSVDofM using:\n";
606   std::cout << "  z: " << col0.transpose() << "\n";
607   std::cout << "  d: " << diag.transpose() << "\n";
608 #endif
609 
610   // Compute singVals, shifts, and mus
611   computeSingVals(col0, diag, perm, singVals, shifts, mus);
612 
613 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
614   std::cout << "  j:        " << (m_computed.block(firstCol, firstCol, n, n)).jacobiSvd().singularValues().transpose().reverse() << "\n\n";
615   std::cout << "  sing-val: " << singVals.transpose() << "\n";
616   std::cout << "  mu:       " << mus.transpose() << "\n";
617   std::cout << "  shift:    " << shifts.transpose() << "\n";
618 
619   {
620     Index actual_n = n;
621     while(actual_n>1 && abs(col0(actual_n-1))<considerZero) --actual_n;
622     std::cout << "\n\n    mus:    " << mus.head(actual_n).transpose() << "\n\n";
623     std::cout << "    check1 (expect0) : " << ((singVals.array()-(shifts+mus)) / singVals.array()).head(actual_n).transpose() << "\n\n";
624     std::cout << "    check2 (>0)      : " << ((singVals.array()-diag) / singVals.array()).head(actual_n).transpose() << "\n\n";
625     std::cout << "    check3 (>0)      : " << ((diag.segment(1,actual_n-1)-singVals.head(actual_n-1).array()) / singVals.head(actual_n-1).array()).transpose() << "\n\n\n";
626     std::cout << "    check4 (>0)      : " << ((singVals.segment(1,actual_n-1)-singVals.head(actual_n-1))).transpose() << "\n\n\n";
627   }
628 #endif
629 
630 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
631   assert(singVals.allFinite());
632   assert(mus.allFinite());
633   assert(shifts.allFinite());
634 #endif
635 
636   // Compute zhat
637   perturbCol0(col0, diag, perm, singVals, shifts, mus, zhat);
638 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
639   std::cout << "  zhat: " << zhat.transpose() << "\n";
640 #endif
641 
642 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
643   assert(zhat.allFinite());
644 #endif
645 
646   computeSingVecs(zhat, diag, perm, singVals, shifts, mus, U, V);
647 
648 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
649   std::cout << "U^T U: " << (U.transpose() * U - MatrixXr(MatrixXr::Identity(U.cols(),U.cols()))).norm() << "\n";
650   std::cout << "V^T V: " << (V.transpose() * V - MatrixXr(MatrixXr::Identity(V.cols(),V.cols()))).norm() << "\n";
651 #endif
652 
653 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
654   assert(U.allFinite());
655   assert(V.allFinite());
656   assert((U.transpose() * U - MatrixXr(MatrixXr::Identity(U.cols(),U.cols()))).norm() < 1e-14 * n);
657   assert((V.transpose() * V - MatrixXr(MatrixXr::Identity(V.cols(),V.cols()))).norm() < 1e-14 * n);
658   assert(m_naiveU.allFinite());
659   assert(m_naiveV.allFinite());
660   assert(m_computed.allFinite());
661 #endif
662 
663   // Because of deflation, the singular values might not be completely sorted.
664   // Fortunately, reordering them is a O(n) problem
665   for(Index i=0; i<actual_n-1; ++i)
666   {
667     if(singVals(i)>singVals(i+1))
668     {
669       using std::swap;
670       swap(singVals(i),singVals(i+1));
671       U.col(i).swap(U.col(i+1));
672       if(m_compV) V.col(i).swap(V.col(i+1));
673     }
674   }
675 
676   // Reverse order so that singular values in increased order
677   // Because of deflation, the zeros singular-values are already at the end
678   singVals.head(actual_n).reverseInPlace();
679   U.leftCols(actual_n).rowwise().reverseInPlace();
680   if (m_compV) V.leftCols(actual_n).rowwise().reverseInPlace();
681 
682 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
683   JacobiSVD<MatrixXr> jsvd(m_computed.block(firstCol, firstCol, n, n) );
684   std::cout << "  * j:        " << jsvd.singularValues().transpose() << "\n\n";
685   std::cout << "  * sing-val: " << singVals.transpose() << "\n";
686 //   std::cout << "  * err:      " << ((jsvd.singularValues()-singVals)>1e-13*singVals.norm()).transpose() << "\n";
687 #endif
688 }
689 
690 template <typename MatrixType>
691 typename BDCSVD<MatrixType>::RealScalar BDCSVD<MatrixType>::secularEq(RealScalar mu, const ArrayRef& col0, const ArrayRef& diag, const IndicesRef &perm, const ArrayRef& diagShifted, RealScalar shift)
692 {
693   Index m = perm.size();
694   RealScalar res = 1;
695   for(Index i=0; i<m; ++i)
696   {
697     Index j = perm(i);
698     res += numext::abs2(col0(j)) / ((diagShifted(j) - mu) * (diag(j) + shift + mu));
699   }
700   return res;
701 
702 }
703 
704 template <typename MatrixType>
705 void BDCSVD<MatrixType>::computeSingVals(const ArrayRef& col0, const ArrayRef& diag, const IndicesRef &perm,
706                                          VectorType& singVals, ArrayRef shifts, ArrayRef mus)
707 {
708   using std::abs;
709   using std::swap;
710 
711   Index n = col0.size();
712   Index actual_n = n;
713   while(actual_n>1 && col0(actual_n-1)==0) --actual_n;
714 
715   for (Index k = 0; k < n; ++k)
716   {
717     if (col0(k) == 0 || actual_n==1)
718     {
719       // if col0(k) == 0, then entry is deflated, so singular value is on diagonal
720       // if actual_n==1, then the deflated problem is already diagonalized
721       singVals(k) = k==0 ? col0(0) : diag(k);
722       mus(k) = 0;
723       shifts(k) = k==0 ? col0(0) : diag(k);
724       continue;
725     }
726 
727     // otherwise, use secular equation to find singular value
728     RealScalar left = diag(k);
729     RealScalar right; // was: = (k != actual_n-1) ? diag(k+1) : (diag(actual_n-1) + col0.matrix().norm());
730     if(k==actual_n-1)
731       right = (diag(actual_n-1) + col0.matrix().norm());
732     else
733     {
734       // Skip deflated singular values
735       Index l = k+1;
736       while(col0(l)==0) { ++l; eigen_internal_assert(l<actual_n); }
737       right = diag(l);
738     }
739 
740     // first decide whether it's closer to the left end or the right end
741     RealScalar mid = left + (right-left) / 2;
742     RealScalar fMid = secularEq(mid, col0, diag, perm, diag, 0);
743 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
744     std::cout << right-left << "\n";
745     std::cout << "fMid = " << fMid << " " << secularEq(mid-left, col0, diag, perm, diag-left, left) << " " << secularEq(mid-right, col0, diag, perm, diag-right, right)   << "\n";
746     std::cout << "     = " << secularEq(0.1*(left+right), col0, diag, perm, diag, 0)
747               << " "       << secularEq(0.2*(left+right), col0, diag, perm, diag, 0)
748               << " "       << secularEq(0.3*(left+right), col0, diag, perm, diag, 0)
749               << " "       << secularEq(0.4*(left+right), col0, diag, perm, diag, 0)
750               << " "       << secularEq(0.49*(left+right), col0, diag, perm, diag, 0)
751               << " "       << secularEq(0.5*(left+right), col0, diag, perm, diag, 0)
752               << " "       << secularEq(0.51*(left+right), col0, diag, perm, diag, 0)
753               << " "       << secularEq(0.6*(left+right), col0, diag, perm, diag, 0)
754               << " "       << secularEq(0.7*(left+right), col0, diag, perm, diag, 0)
755               << " "       << secularEq(0.8*(left+right), col0, diag, perm, diag, 0)
756               << " "       << secularEq(0.9*(left+right), col0, diag, perm, diag, 0) << "\n";
757 #endif
758     RealScalar shift = (k == actual_n-1 || fMid > 0) ? left : right;
759 
760     // measure everything relative to shift
761     Map<ArrayXr> diagShifted(m_workspace.data()+4*n, n);
762     diagShifted = diag - shift;
763 
764     // initial guess
765     RealScalar muPrev, muCur;
766     if (shift == left)
767     {
768       muPrev = (right - left) * RealScalar(0.1);
769       if (k == actual_n-1) muCur = right - left;
770       else                 muCur = (right - left) * RealScalar(0.5);
771     }
772     else
773     {
774       muPrev = -(right - left) * RealScalar(0.1);
775       muCur = -(right - left) * RealScalar(0.5);
776     }
777 
778     RealScalar fPrev = secularEq(muPrev, col0, diag, perm, diagShifted, shift);
779     RealScalar fCur = secularEq(muCur, col0, diag, perm, diagShifted, shift);
780     if (abs(fPrev) < abs(fCur))
781     {
782       swap(fPrev, fCur);
783       swap(muPrev, muCur);
784     }
785 
786     // rational interpolation: fit a function of the form a / mu + b through the two previous
787     // iterates and use its zero to compute the next iterate
788     bool useBisection = fPrev*fCur>0;
789     while (fCur!=0 && abs(muCur - muPrev) > 8 * NumTraits<RealScalar>::epsilon() * numext::maxi<RealScalar>(abs(muCur), abs(muPrev)) && abs(fCur - fPrev)>NumTraits<RealScalar>::epsilon() && !useBisection)
790     {
791       ++m_numIters;
792 
793       // Find a and b such that the function f(mu) = a / mu + b matches the current and previous samples.
794       RealScalar a = (fCur - fPrev) / (1/muCur - 1/muPrev);
795       RealScalar b = fCur - a / muCur;
796       // And find mu such that f(mu)==0:
797       RealScalar muZero = -a/b;
798       RealScalar fZero = secularEq(muZero, col0, diag, perm, diagShifted, shift);
799 
800       muPrev = muCur;
801       fPrev = fCur;
802       muCur = muZero;
803       fCur = fZero;
804 
805 
806       if (shift == left  && (muCur < 0 || muCur > right - left)) useBisection = true;
807       if (shift == right && (muCur < -(right - left) || muCur > 0)) useBisection = true;
808       if (abs(fCur)>abs(fPrev)) useBisection = true;
809     }
810 
811     // fall back on bisection method if rational interpolation did not work
812     if (useBisection)
813     {
814 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
815       std::cout << "useBisection for k = " << k << ", actual_n = " << actual_n << "\n";
816 #endif
817       RealScalar leftShifted, rightShifted;
818       if (shift == left)
819       {
820         leftShifted = (std::numeric_limits<RealScalar>::min)();
821         // I don't understand why the case k==0 would be special there:
822         // if (k == 0) rightShifted = right - left; else
823         rightShifted = (k==actual_n-1) ? right : ((right - left) * RealScalar(0.6)); // theoretically we can take 0.5, but let's be safe
824       }
825       else
826       {
827         leftShifted = -(right - left) * RealScalar(0.6);
828         rightShifted = -(std::numeric_limits<RealScalar>::min)();
829       }
830 
831       RealScalar fLeft = secularEq(leftShifted, col0, diag, perm, diagShifted, shift);
832 
833 #if defined EIGEN_INTERNAL_DEBUGGING || defined EIGEN_BDCSVD_DEBUG_VERBOSE
834       RealScalar fRight = secularEq(rightShifted, col0, diag, perm, diagShifted, shift);
835 #endif
836 
837 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
838       if(!(fLeft * fRight<0))
839       {
840         std::cout << "fLeft: " << leftShifted << " - " << diagShifted.head(10).transpose()  << "\n ; " << bool(left==shift) << " " << (left-shift) << "\n";
841         std::cout << k << " : " <<  fLeft << " * " << fRight << " == " << fLeft * fRight << "  ;  " << left << " - " << right << " -> " <<  leftShifted << " " << rightShifted << "   shift=" << shift << "\n";
842       }
843 #endif
844       eigen_internal_assert(fLeft * fRight < 0);
845 
846       while (rightShifted - leftShifted > 2 * NumTraits<RealScalar>::epsilon() * numext::maxi<RealScalar>(abs(leftShifted), abs(rightShifted)))
847       {
848         RealScalar midShifted = (leftShifted + rightShifted) / 2;
849         fMid = secularEq(midShifted, col0, diag, perm, diagShifted, shift);
850         if (fLeft * fMid < 0)
851         {
852           rightShifted = midShifted;
853         }
854         else
855         {
856           leftShifted = midShifted;
857           fLeft = fMid;
858         }
859       }
860 
861       muCur = (leftShifted + rightShifted) / 2;
862     }
863 
864     singVals[k] = shift + muCur;
865     shifts[k] = shift;
866     mus[k] = muCur;
867 
868     // perturb singular value slightly if it equals diagonal entry to avoid division by zero later
869     // (deflation is supposed to avoid this from happening)
870     // - this does no seem to be necessary anymore -
871 //     if (singVals[k] == left) singVals[k] *= 1 + NumTraits<RealScalar>::epsilon();
872 //     if (singVals[k] == right) singVals[k] *= 1 - NumTraits<RealScalar>::epsilon();
873   }
874 }
875 
876 
877 // zhat is perturbation of col0 for which singular vectors can be computed stably (see Section 3.1)
878 template <typename MatrixType>
879 void BDCSVD<MatrixType>::perturbCol0
880    (const ArrayRef& col0, const ArrayRef& diag, const IndicesRef &perm, const VectorType& singVals,
881     const ArrayRef& shifts, const ArrayRef& mus, ArrayRef zhat)
882 {
883   using std::sqrt;
884   Index n = col0.size();
885   Index m = perm.size();
886   if(m==0)
887   {
888     zhat.setZero();
889     return;
890   }
891   Index last = perm(m-1);
892   // The offset permits to skip deflated entries while computing zhat
893   for (Index k = 0; k < n; ++k)
894   {
895     if (col0(k) == 0) // deflated
896       zhat(k) = 0;
897     else
898     {
899       // see equation (3.6)
900       RealScalar dk = diag(k);
901       RealScalar prod = (singVals(last) + dk) * (mus(last) + (shifts(last) - dk));
902 
903       for(Index l = 0; l<m; ++l)
904       {
905         Index i = perm(l);
906         if(i!=k)
907         {
908           Index j = i<k ? i : perm(l-1);
909           prod *= ((singVals(j)+dk) / ((diag(i)+dk))) * ((mus(j)+(shifts(j)-dk)) / ((diag(i)-dk)));
910 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
911           if(i!=k && std::abs(((singVals(j)+dk)*(mus(j)+(shifts(j)-dk)))/((diag(i)+dk)*(diag(i)-dk)) - 1) > 0.9 )
912             std::cout << "     " << ((singVals(j)+dk)*(mus(j)+(shifts(j)-dk)))/((diag(i)+dk)*(diag(i)-dk)) << " == (" << (singVals(j)+dk) << " * " << (mus(j)+(shifts(j)-dk))
913                        << ") / (" << (diag(i)+dk) << " * " << (diag(i)-dk) << ")\n";
914 #endif
915         }
916       }
917 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
918       std::cout << "zhat(" << k << ") =  sqrt( " << prod << ")  ;  " << (singVals(last) + dk) << " * " << mus(last) + shifts(last) << " - " << dk << "\n";
919 #endif
920       RealScalar tmp = sqrt(prod);
921       zhat(k) = col0(k) > 0 ? tmp : -tmp;
922     }
923   }
924 }
925 
926 // compute singular vectors
927 template <typename MatrixType>
928 void BDCSVD<MatrixType>::computeSingVecs
929    (const ArrayRef& zhat, const ArrayRef& diag, const IndicesRef &perm, const VectorType& singVals,
930     const ArrayRef& shifts, const ArrayRef& mus, MatrixXr& U, MatrixXr& V)
931 {
932   Index n = zhat.size();
933   Index m = perm.size();
934 
935   for (Index k = 0; k < n; ++k)
936   {
937     if (zhat(k) == 0)
938     {
939       U.col(k) = VectorType::Unit(n+1, k);
940       if (m_compV) V.col(k) = VectorType::Unit(n, k);
941     }
942     else
943     {
944       U.col(k).setZero();
945       for(Index l=0;l<m;++l)
946       {
947         Index i = perm(l);
948         U(i,k) = zhat(i)/(((diag(i) - shifts(k)) - mus(k)) )/( (diag(i) + singVals[k]));
949       }
950       U(n,k) = 0;
951       U.col(k).normalize();
952 
953       if (m_compV)
954       {
955         V.col(k).setZero();
956         for(Index l=1;l<m;++l)
957         {
958           Index i = perm(l);
959           V(i,k) = diag(i) * zhat(i) / (((diag(i) - shifts(k)) - mus(k)) )/( (diag(i) + singVals[k]));
960         }
961         V(0,k) = -1;
962         V.col(k).normalize();
963       }
964     }
965   }
966   U.col(n) = VectorType::Unit(n+1, n);
967 }
968 
969 
970 // page 12_13
971 // i >= 1, di almost null and zi non null.
972 // We use a rotation to zero out zi applied to the left of M
973 template <typename MatrixType>
974 void BDCSVD<MatrixType>::deflation43(Index firstCol, Index shift, Index i, Index size)
975 {
976   using std::abs;
977   using std::sqrt;
978   using std::pow;
979   Index start = firstCol + shift;
980   RealScalar c = m_computed(start, start);
981   RealScalar s = m_computed(start+i, start);
982   RealScalar r = sqrt(numext::abs2(c) + numext::abs2(s));
983   if (r == 0)
984   {
985     m_computed(start+i, start+i) = 0;
986     return;
987   }
988   m_computed(start,start) = r;
989   m_computed(start+i, start) = 0;
990   m_computed(start+i, start+i) = 0;
991 
992   JacobiRotation<RealScalar> J(c/r,-s/r);
993   if (m_compU)  m_naiveU.middleRows(firstCol, size+1).applyOnTheRight(firstCol, firstCol+i, J);
994   else          m_naiveU.applyOnTheRight(firstCol, firstCol+i, J);
995 }// end deflation 43
996 
997 
998 // page 13
999 // i,j >= 1, i!=j and |di - dj| < epsilon * norm2(M)
1000 // We apply two rotations to have zj = 0;
1001 // TODO deflation44 is still broken and not properly tested
1002 template <typename MatrixType>
1003 void BDCSVD<MatrixType>::deflation44(Index firstColu , Index firstColm, Index firstRowW, Index firstColW, Index i, Index j, Index size)
1004 {
1005   using std::abs;
1006   using std::sqrt;
1007   using std::conj;
1008   using std::pow;
1009   RealScalar c = m_computed(firstColm+i, firstColm);
1010   RealScalar s = m_computed(firstColm+j, firstColm);
1011   RealScalar r = sqrt(numext::abs2(c) + numext::abs2(s));
1012 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
1013   std::cout << "deflation 4.4: " << i << "," << j << " -> " << c << " " << s << " " << r << " ; "
1014     << m_computed(firstColm + i-1, firstColm)  << " "
1015     << m_computed(firstColm + i, firstColm)  << " "
1016     << m_computed(firstColm + i+1, firstColm) << " "
1017     << m_computed(firstColm + i+2, firstColm) << "\n";
1018   std::cout << m_computed(firstColm + i-1, firstColm + i-1)  << " "
1019     << m_computed(firstColm + i, firstColm+i)  << " "
1020     << m_computed(firstColm + i+1, firstColm+i+1) << " "
1021     << m_computed(firstColm + i+2, firstColm+i+2) << "\n";
1022 #endif
1023   if (r==0)
1024   {
1025     m_computed(firstColm + i, firstColm + i) = m_computed(firstColm + j, firstColm + j);
1026     return;
1027   }
1028   c/=r;
1029   s/=r;
1030   m_computed(firstColm + i, firstColm) = r;
1031   m_computed(firstColm + j, firstColm + j) = m_computed(firstColm + i, firstColm + i);
1032   m_computed(firstColm + j, firstColm) = 0;
1033 
1034   JacobiRotation<RealScalar> J(c,-s);
1035   if (m_compU)  m_naiveU.middleRows(firstColu, size+1).applyOnTheRight(firstColu + i, firstColu + j, J);
1036   else          m_naiveU.applyOnTheRight(firstColu+i, firstColu+j, J);
1037   if (m_compV)  m_naiveV.middleRows(firstRowW, size).applyOnTheRight(firstColW + i, firstColW + j, J);
1038 }// end deflation 44
1039 
1040 
1041 // acts on block from (firstCol+shift, firstCol+shift) to (lastCol+shift, lastCol+shift) [inclusive]
1042 template <typename MatrixType>
1043 void BDCSVD<MatrixType>::deflation(Index firstCol, Index lastCol, Index k, Index firstRowW, Index firstColW, Index shift)
1044 {
1045   using std::sqrt;
1046   using std::abs;
1047   const Index length = lastCol + 1 - firstCol;
1048 
1049   Block<MatrixXr,Dynamic,1> col0(m_computed, firstCol+shift, firstCol+shift, length, 1);
1050   Diagonal<MatrixXr> fulldiag(m_computed);
1051   VectorBlock<Diagonal<MatrixXr>,Dynamic> diag(fulldiag, firstCol+shift, length);
1052 
1053   const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
1054   RealScalar maxDiag = diag.tail((std::max)(Index(1),length-1)).cwiseAbs().maxCoeff();
1055   RealScalar epsilon_strict = numext::maxi<RealScalar>(considerZero,NumTraits<RealScalar>::epsilon() * maxDiag);
1056   RealScalar epsilon_coarse = 8 * NumTraits<RealScalar>::epsilon() * numext::maxi<RealScalar>(col0.cwiseAbs().maxCoeff(), maxDiag);
1057 
1058 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1059   assert(m_naiveU.allFinite());
1060   assert(m_naiveV.allFinite());
1061   assert(m_computed.allFinite());
1062 #endif
1063 
1064 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
1065   std::cout << "\ndeflate:" << diag.head(k+1).transpose() << "  |  " << diag.segment(k+1,length-k-1).transpose() << "\n";
1066 #endif
1067 
1068   //condition 4.1
1069   if (diag(0) < epsilon_coarse)
1070   {
1071 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
1072     std::cout << "deflation 4.1, because " << diag(0) << " < " << epsilon_coarse << "\n";
1073 #endif
1074     diag(0) = epsilon_coarse;
1075   }
1076 
1077   //condition 4.2
1078   for (Index i=1;i<length;++i)
1079     if (abs(col0(i)) < epsilon_strict)
1080     {
1081 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
1082       std::cout << "deflation 4.2, set z(" << i << ") to zero because " << abs(col0(i)) << " < " << epsilon_strict << "  (diag(" << i << ")=" << diag(i) << ")\n";
1083 #endif
1084       col0(i) = 0;
1085     }
1086 
1087   //condition 4.3
1088   for (Index i=1;i<length; i++)
1089     if (diag(i) < epsilon_coarse)
1090     {
1091 #ifdef  EIGEN_BDCSVD_DEBUG_VERBOSE
1092       std::cout << "deflation 4.3, cancel z(" << i << ")=" << col0(i) << " because diag(" << i << ")=" << diag(i) << " < " << epsilon_coarse << "\n";
1093 #endif
1094       deflation43(firstCol, shift, i, length);
1095     }
1096 
1097 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1098   assert(m_naiveU.allFinite());
1099   assert(m_naiveV.allFinite());
1100   assert(m_computed.allFinite());
1101 #endif
1102 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1103   std::cout << "to be sorted: " << diag.transpose() << "\n\n";
1104 #endif
1105   {
1106     // Check for total deflation
1107     // If we have a total deflation, then we have to consider col0(0)==diag(0) as a singular value during sorting
1108     bool total_deflation = (col0.tail(length-1).array()<considerZero).all();
1109 
1110     // Sort the diagonal entries, since diag(1:k-1) and diag(k:length) are already sorted, let's do a sorted merge.
1111     // First, compute the respective permutation.
1112     Index *permutation = m_workspaceI.data();
1113     {
1114       permutation[0] = 0;
1115       Index p = 1;
1116 
1117       // Move deflated diagonal entries at the end.
1118       for(Index i=1; i<length; ++i)
1119         if(abs(diag(i))<considerZero)
1120           permutation[p++] = i;
1121 
1122       Index i=1, j=k+1;
1123       for( ; p < length; ++p)
1124       {
1125              if (i > k)             permutation[p] = j++;
1126         else if (j >= length)       permutation[p] = i++;
1127         else if (diag(i) < diag(j)) permutation[p] = j++;
1128         else                        permutation[p] = i++;
1129       }
1130     }
1131 
1132     // If we have a total deflation, then we have to insert diag(0) at the right place
1133     if(total_deflation)
1134     {
1135       for(Index i=1; i<length; ++i)
1136       {
1137         Index pi = permutation[i];
1138         if(abs(diag(pi))<considerZero || diag(0)<diag(pi))
1139           permutation[i-1] = permutation[i];
1140         else
1141         {
1142           permutation[i-1] = 0;
1143           break;
1144         }
1145       }
1146     }
1147 
1148     // Current index of each col, and current column of each index
1149     Index *realInd = m_workspaceI.data()+length;
1150     Index *realCol = m_workspaceI.data()+2*length;
1151 
1152     for(int pos = 0; pos< length; pos++)
1153     {
1154       realCol[pos] = pos;
1155       realInd[pos] = pos;
1156     }
1157 
1158     for(Index i = total_deflation?0:1; i < length; i++)
1159     {
1160       const Index pi = permutation[length - (total_deflation ? i+1 : i)];
1161       const Index J = realCol[pi];
1162 
1163       using std::swap;
1164       // swap diagonal and first column entries:
1165       swap(diag(i), diag(J));
1166       if(i!=0 && J!=0) swap(col0(i), col0(J));
1167 
1168       // change columns
1169       if (m_compU) m_naiveU.col(firstCol+i).segment(firstCol, length + 1).swap(m_naiveU.col(firstCol+J).segment(firstCol, length + 1));
1170       else         m_naiveU.col(firstCol+i).segment(0, 2)                .swap(m_naiveU.col(firstCol+J).segment(0, 2));
1171       if (m_compV) m_naiveV.col(firstColW + i).segment(firstRowW, length).swap(m_naiveV.col(firstColW + J).segment(firstRowW, length));
1172 
1173       //update real pos
1174       const Index realI = realInd[i];
1175       realCol[realI] = J;
1176       realCol[pi] = i;
1177       realInd[J] = realI;
1178       realInd[i] = pi;
1179     }
1180   }
1181 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1182   std::cout << "sorted: " << diag.transpose().format(bdcsvdfmt) << "\n";
1183   std::cout << "      : " << col0.transpose() << "\n\n";
1184 #endif
1185 
1186   //condition 4.4
1187   {
1188     Index i = length-1;
1189     while(i>0 && (abs(diag(i))<considerZero || abs(col0(i))<considerZero)) --i;
1190     for(; i>1;--i)
1191        if( (diag(i) - diag(i-1)) < NumTraits<RealScalar>::epsilon()*maxDiag )
1192       {
1193 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1194         std::cout << "deflation 4.4 with i = " << i << " because " << (diag(i) - diag(i-1)) << " < " << NumTraits<RealScalar>::epsilon()*diag(i) << "\n";
1195 #endif
1196         eigen_internal_assert(abs(diag(i) - diag(i-1))<epsilon_coarse && " diagonal entries are not properly sorted");
1197         deflation44(firstCol, firstCol + shift, firstRowW, firstColW, i-1, i, length);
1198       }
1199   }
1200 
1201 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1202   for(Index j=2;j<length;++j)
1203     assert(diag(j-1)<=diag(j) || abs(diag(j))<considerZero);
1204 #endif
1205 
1206 #ifdef EIGEN_BDCSVD_SANITY_CHECKS
1207   assert(m_naiveU.allFinite());
1208   assert(m_naiveV.allFinite());
1209   assert(m_computed.allFinite());
1210 #endif
1211 }//end deflation
1212 
1213 #ifndef __CUDACC__
1214 /** \svd_module
1215   *
1216   * \return the singular value decomposition of \c *this computed by Divide & Conquer algorithm
1217   *
1218   * \sa class BDCSVD
1219   */
1220 template<typename Derived>
1221 BDCSVD<typename MatrixBase<Derived>::PlainObject>
1222 MatrixBase<Derived>::bdcSvd(unsigned int computationOptions) const
1223 {
1224   return BDCSVD<PlainObject>(*this, computationOptions);
1225 }
1226 #endif
1227 
1228 } // end namespace Eigen
1229 
1230 #endif
1231