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