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