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