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